source: LMDZ6/branches/Amaury_dev/libf/dynphy_lonlat/phylmd/limit_netcdf.f90 @ 5218

Last change on this file since 5218 was 5159, checked in by abarral, 3 months ago

Put dimensions.h and paramet.h into modules

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 39.1 KB
Line 
1MODULE limit
2
3  !*******************************************************************************
4  ! Author : L. Fairhead, 27/01/94
5  !-------------------------------------------------------------------------------
6  ! Purpose: Boundary conditions files building for new model using climatologies.
7  !          Both grids have to be regular.
8  !-------------------------------------------------------------------------------
9  ! Note: This routine is designed to work for Earth
10  !-------------------------------------------------------------------------------
11  ! Modification history:
12  !  * 23/03/1994: Z. X. Li
13  !  *    09/1999: L. Fairhead (netcdf reading in LMDZ.3.3)
14  !  *    07/2001: P. Le Van
15  !  *    11/2009: L. Guez     (ozone day & night climatos, see etat0_netcdf.F90)
16  !  *    12/2009: D. Cugnet   (f77->f90, calendars, files from coupled runs)
17  !-------------------------------------------------------------------------------
18
19  USE ioipsl, ONLY: flininfo, flinopen, flinget, flinclo
20  USE lmdz_assert_eq, ONLY: assert_eq
21  USE cal_tools_m, ONLY: year_len, mid_month
22  USE conf_dat_m, ONLY: conf_dat2d, conf_dat3d
23  USE dimphy, ONLY: klon, zmasq
24  USE lmdz_geometry, ONLY: longitude_deg, latitude_deg
25  USE phys_state_var_mod, ONLY: pctsrf
26  USE control_mod, ONLY: anneeref
27  USE init_ssrf_m, ONLY: start_init_subsurf
28
29  INTEGER, PARAMETER :: ns = 256
30  CHARACTER(LEN = ns), PARAMETER :: &
31          fsst(5) = ['amipbc_sst_1x1.nc   ', 'amip_sst_1x1.nc     ', 'cpl_atm_sst.nc      '&
32                  , 'histmth_sst.nc      ', 'sstk.nc             '], &
33          fsic(5) = ['amipbc_sic_1x1.nc   ', 'amip_sic_1x1.nc     ', 'cpl_atm_sic.nc      '&
34                  , 'histmth_sic.nc      ', 'ci.nc               '], &
35          vsst(5) = ['tosbcs    ', 'tos       ', 'SISUTESW  ', 'tsol_oce  ', 'sstk      '], &
36          vsic(5) = ['sicbcs    ', 'sic       ', 'SIICECOV  ', 'pourc_sic ', 'ci        '], &
37          frugo = 'Rugos.nc  ', falbe = 'Albedo.nc ', frelf = 'Relief.nc ', &
38          vrug = 'RUGOS     ', valb = 'ALBEDO    ', vrel = 'RELIEF    ', &
39          DegK(11) = ['degK          ', 'degree_K      ', 'degreeK       ', 'deg_K         '&
40                  , 'degsK         ', 'degrees_K     ', 'degreesK      ', 'degs_K        '&
41                  , 'degree_kelvin ', 'degrees_kelvin', 'K             '], &
42          DegC(10) = ['degC          ', 'degree_C      ', 'degreeC       ', 'deg_C         '&
43                  , 'degsC         ', 'degrees_C     ', 'degreesC      ', 'degs_C        '&
44                  , 'degree_Celsius', 'celsius       '], &
45          Perc(2) = ['%             ', 'percent       '], &
46          Frac(2) = ['1.0           ', '1             ']
47
48CONTAINS
49
50  !-------------------------------------------------------------------------------
51
52  SUBROUTINE limit_netcdf(masque, phis, extrap)
53
54    !-------------------------------------------------------------------------------
55    ! Author : L. Fairhead, 27/01/94
56    !-------------------------------------------------------------------------------
57    ! Purpose: Boundary conditions files building for new model using climatologies.
58    !          Both grids have to be regular.
59    !-------------------------------------------------------------------------------
60    ! Note: This routine is designed to work for Earth
61    !-------------------------------------------------------------------------------
62    ! Modification history:
63    !  * 23/03/1994: Z. X. Li
64    !  *    09/1999: L. Fairhead (netcdf reading in LMDZ.3.3)
65    !  *    07/2001: P. Le Van
66    !  *    11/2009: L. Guez     (ozone day & night climatos, see etat0_netcdf.F90)
67    !  *    12/2009: D. Cugnet   (f77->f90, calendars, files from coupled runs)
68    !  *    04/2016: D. Cugnet   (12/14 recs SST/SIC files: cyclic/interannual runs)
69    !  *    05/2017: D. Cugnet   (linear time interpolation for BCS files)
70    !-------------------------------------------------------------------------------
71    USE indice_sol_mod
72    USE netcdf, ONLY: nf90_open, nf90_create, nf90_close, &
73            nf90_def_dim, nf90_def_var, nf90_put_var, nf90_put_att, &
74            nf90_noerr, nf90_nowrite, nf90_global, &
75            nf90_clobber, nf90_enddef, nf90_unlimited, nf90_float, &
76            nf90_64bit_offset
77    USE lmdz_cppkeys_wrapper, ONLY: nf90_format
78    USE inter_barxy_m, ONLY: inter_barxy
79    USE netcdf95, ONLY: nf95_def_var, nf95_put_att, nf95_put_var
80    USE comconst_mod, ONLY: pi
81    USE phys_cal_mod, ONLY: calend
82    USE lmdz_abort_physic, ONLY: abort_physic
83    USE lmdz_iniprint, ONLY: lunout, prt_level
84    USE lmdz_comgeom2
85
86  USE lmdz_dimensions, ONLY: iim, jjm, llm, ndm
87  USE lmdz_paramet
88    IMPLICIT NONE
89    !-------------------------------------------------------------------------------
90    ! Arguments:
91
92
93    REAL, DIMENSION(iip1, jjp1), INTENT(INOUT) :: masque ! land mask
94    REAL, DIMENSION(iip1, jjp1), INTENT(INOUT) :: phis   ! ground geopotential
95    LOGICAL, INTENT(IN) :: extrap ! SST extrapolation flag
96    !-------------------------------------------------------------------------------
97    ! Local variables:
98
99    !--- INPUT NETCDF FILES AND VARIABLES NAMES ------------------------------------
100    CHARACTER(LEN = ns) :: icefile, sstfile, fnam, varname
101
102    !--- OUTPUT VARIABLES FOR NETCDF FILE ------------------------------------------
103    REAL :: fi_ice(klon)
104    REAL, POINTER :: phy_rug(:, :) => NULL(), phy_ice(:, :) => NULL()
105    REAL, POINTER :: phy_sst(:, :) => NULL(), phy_alb(:, :) => NULL()
106    REAL, ALLOCATABLE :: phy_bil(:, :), pctsrf_t(:, :, :)
107    INTEGER :: nbad
108
109    !--- VARIABLES FOR OUTPUT FILE WRITING -----------------------------------------
110    INTEGER :: nid, ndim, ntim, k, dims(2), ix_sic, ix_sst
111    INTEGER :: id_tim, id_SST, id_BILS, id_RUG, id_ALB
112    INTEGER :: id_FOCE, id_FSIC, id_FTER, id_FLIC, varid_longitude, varid_latitude
113    INTEGER :: ndays                   !--- Depending on the output calendar
114    CHARACTER(LEN = ns) :: str
115
116    !--- INITIALIZATIONS -----------------------------------------------------------
117    CALL inigeom
118
119    !--- MASK, GROUND GEOPOT. & SUBSURFACES COMPUTATION (IN CASE ok_etat0==.FALSE.)
120    IF(ALL(masque==-99999.)) THEN
121      CALL start_init_orog0(rlonv, rlatu, phis, masque)
122      CALL gr_dyn_fi(1, iip1, jjp1, klon, masque, zmasq)          !--- To physical grid
123      ALLOCATE(pctsrf(klon, nbsrf))
124      CALL start_init_subsurf(.FALSE.)
125      !--- TO MATCH EXACTLY WHAT WOULD BE DONE IN etat0phys_netcdf
126      WHERE(masque(:, :)<EPSFRA) masque(:, :) = 0.
127      WHERE(1. - masque(:, :)<EPSFRA) masque(:, :) = 1.
128    END IF
129
130    !--- Beware: anneeref (from gcm.def) is used to determine output time sampling
131    ndays = year_len(anneeref)
132
133    !--- RUGOSITY TREATMENT --------------------------------------------------------
134    CALL msg(0, ""); CALL msg(0, " *** TRAITEMENT DE LA RUGOSITE ***")
135    CALL get_2Dfield(frugo, vrug, 'RUG', ndays, phy_rug, mask = masque(1:iim, :))
136
137    !--- OCEAN TREATMENT -----------------------------------------------------------
138    CALL msg(0, ""); CALL msg(0, " *** TRAITEMENT DE LA GLACE OCEANIQUE ***")
139
140    ! Input SIC file selection
141    ! Open file only to test if available
142    DO ix_sic = 1, SIZE(fsic)
143      IF (nf90_open(TRIM(fsic(ix_sic)), nf90_nowrite, nid)==nf90_noerr) THEN
144        icefile = fsic(ix_sic); varname = vsic(ix_sic); EXIT
145      END IF
146    END DO
147    IF(ix_sic==SIZE(fsic) + 1) THEN
148      WRITE(lunout, *) 'ERROR! No sea-ice input file was found.'
149      WRITE(lunout, *) 'One of following files must be available : '
150      DO k = 1, SIZE(fsic); WRITE(lunout, *) TRIM(fsic(k));
151      END DO
152      CALL abort_physic('limit_netcdf', 'No sea-ice file was found', 1)
153    END IF
154    CALL ncerr(nf90_close(nid), icefile)
155    CALL msg(0, 'Fichier choisi pour la glace de mer:' // TRIM(icefile))
156
157    CALL get_2Dfield(icefile, varname, 'SIC', ndays, phy_ice)
158
159    ALLOCATE(pctsrf_t(klon, nbsrf, ndays))
160    DO k = 1, ndays
161      fi_ice = phy_ice(:, k)
162      WHERE(fi_ice>=1.0) fi_ice = 1.0
163      WHERE(fi_ice<EPSFRA) fi_ice = 0.0
164      pctsrf_t(:, is_ter, k) = pctsrf(:, is_ter)       ! land soil
165      pctsrf_t(:, is_lic, k) = pctsrf(:, is_lic)       ! land ice
166      SELECT CASE(ix_sic)
167      CASE(3)                                   ! SIC=pICE*(1-LIC-TER) (CPL)
168        pctsrf_t(:, is_sic, k) = fi_ice(:) * (1. - pctsrf(:, is_lic) - pctsrf(:, is_ter))
169      CASE(4)                                   ! SIC=pICE            (HIST)
170        pctsrf_t(:, is_sic, k) = fi_ice(:)
171      CASE DEFAULT                              ! SIC=pICE-LIC   (AMIP,ERAI)
172        pctsrf_t(:, is_sic, k) = fi_ice - pctsrf_t(:, is_lic, k)
173      END SELECT
174      WHERE(pctsrf_t(:, is_sic, k)<=0) pctsrf_t(:, is_sic, k) = 0.
175      WHERE(1.0 - zmasq<EPSFRA)
176        pctsrf_t(:, is_sic, k) = 0.0
177        pctsrf_t(:, is_oce, k) = 0.0
178      ELSEWHERE
179        WHERE(pctsrf_t(:, is_sic, k)>=1.0 - zmasq)
180          pctsrf_t(:, is_sic, k) = 1.0 - zmasq
181          pctsrf_t(:, is_oce, k) = 0.0
182        ELSEWHERE
183          pctsrf_t(:, is_oce, k) = 1.0 - zmasq - pctsrf_t(:, is_sic, k)
184          WHERE(pctsrf_t(:, is_oce, k)<EPSFRA)
185            pctsrf_t(:, is_oce, k) = 0.0
186            pctsrf_t(:, is_sic, k) = 1.0 - zmasq
187          END WHERE
188        END WHERE
189      END WHERE
190      nbad = COUNT(pctsrf_t(:, is_oce, k)<0.0)
191      IF(nbad>0) WRITE(lunout, *) 'pb sous maille pour nb points = ', nbad
192      nbad = COUNT(ABS(SUM(pctsrf_t(:, :, k), DIM = 2) - 1.0)>EPSFRA)
193      IF(nbad>0) WRITE(lunout, *) 'pb sous surface pour nb points = ', nbad
194    END DO
195    DEALLOCATE(phy_ice)
196
197    !--- SST TREATMENT -------------------------------------------------------------
198    CALL msg(0, ""); CALL msg(0, " *** TRAITEMENT DE LA SST ***")
199
200    ! Input SST file selection
201    ! Open file only to test if available
202    DO ix_sst = 1, SIZE(fsst)
203      IF (nf90_open(TRIM(fsst(ix_sst)), nf90_nowrite, nid)==nf90_noerr) THEN
204        sstfile = fsst(ix_sst); varname = vsst(ix_sst); EXIT
205      END IF
206    END DO
207    IF(ix_sst==SIZE(fsst) + 1) THEN
208      WRITE(lunout, *) 'ERROR! No sst input file was found.'
209      WRITE(lunout, *) 'One of following files must be available : '
210      DO k = 1, SIZE(fsst); WRITE(lunout, *) TRIM(fsst(k));
211      END DO
212      CALL abort_physic('limit_netcdf', 'No sst file was found', 1)
213    END IF
214    CALL ncerr(nf90_close(nid), sstfile)
215    CALL msg(0, 'Fichier choisi pour la temperature de mer: ' // TRIM(sstfile))
216
217    CALL get_2Dfield(sstfile, varname, 'SST', ndays, phy_sst, flag = extrap)
218
219    !--- ALBEDO TREATMENT ----------------------------------------------------------
220    CALL msg(0, ""); CALL msg(0, " *** TRAITEMENT DE L'ALBEDO ***")
221    CALL get_2Dfield(falbe, valb, 'ALB', ndays, phy_alb)
222
223    !--- REFERENCE GROUND HEAT FLUX TREATMENT --------------------------------------
224    ALLOCATE(phy_bil(klon, ndays)); phy_bil = 0.0
225
226    !--- OUTPUT FILE WRITING -------------------------------------------------------
227    CALL msg(0, ""); CALL msg(0, ' *** Ecriture du fichier limit : debut ***')
228    fnam = "limit.nc"
229
230    !--- File creation
231    CALL ncerr(nf90_create(fnam, IOR(nf90_clobber, nf90_64bit_offset), nid), fnam)
232    CALL ncerr(nf90_put_att(nid, nf90_global, "title", "Fichier conditions aux limites"), fnam)
233    str = 'File produced using ce0l executable.'
234    str = TRIM(str) // NEW_LINE(' ') // 'Sea Ice Concentration built from'
235    SELECT CASE(ix_sic)
236    CASE(1); str = TRIM(str) // ' Amip mid-month boundary condition (BCS).'
237    CASE(2); str = TRIM(str) // ' Amip monthly mean observations.'
238    CASE(3); str = TRIM(str) // ' IPSL coupled model outputs.'
239    CASE(4); str = TRIM(str) // ' LMDZ model outputs.'
240    CASE(5); str = TRIM(str) // ' ci.nc file.'
241    END SELECT
242    str = TRIM(str) // NEW_LINE(' ') // 'Sea Surface Temperature built from'
243    SELECT CASE(ix_sst)
244    CASE(1); str = TRIM(str) // ' Amip mid-month boundary condition (BCS).'
245    CASE(2); str = TRIM(str) // ' Amip monthly mean observations.'
246    CASE(3); str = TRIM(str) // ' IPSL coupled model outputs.'
247    CASE(4); str = TRIM(str) // ' LMDZ model outputs.'
248    CASE(5); str = TRIM(str) // ' sstk.nc file.'
249    END SELECT
250    CALL ncerr(nf90_put_att(nid, nf90_global, "history", TRIM(str)), fnam)
251
252    !--- Dimensions creation
253    CALL ncerr(nf90_def_dim(nid, "points_physiques", klon, ndim), fnam)
254    CALL ncerr(nf90_def_dim(nid, "time", nf90_unlimited, ntim), fnam)
255
256    dims = [ndim, ntim]
257
258    !--- Variables creation
259    CALL ncerr(nf90_def_var(nid, "TEMPS", nf90_format, [ntim], id_tim), fnam)
260    CALL ncerr(nf90_def_var(nid, "FOCE", nf90_format, dims, id_FOCE), fnam)
261    CALL ncerr(nf90_def_var(nid, "FSIC", nf90_format, dims, id_FSIC), fnam)
262    CALL ncerr(nf90_def_var(nid, "FTER", nf90_format, dims, id_FTER), fnam)
263    CALL ncerr(nf90_def_var(nid, "FLIC", nf90_format, dims, id_FLIC), fnam)
264    CALL ncerr(nf90_def_var(nid, "SST", nf90_format, dims, id_SST), fnam)
265    CALL ncerr(nf90_def_var(nid, "BILS", nf90_format, dims, id_BILS), fnam)
266    CALL ncerr(nf90_def_var(nid, "ALB", nf90_format, dims, id_ALB), fnam)
267    CALL ncerr(nf90_def_var(nid, "RUG", nf90_format, dims, id_RUG), fnam)
268    CALL nf95_def_var(nid, "longitude", nf90_float, ndim, varid_longitude)
269    CALL nf95_def_var(nid, "latitude", nf90_float, ndim, varid_latitude)
270
271    !--- Attributes creation
272    CALL ncerr(nf90_put_att(nid, id_tim, "title", "Jour dans l annee"), fnam)
273    CALL ncerr(nf90_put_att(nid, id_tim, "calendar", calend), fnam)
274    CALL ncerr(nf90_put_att(nid, id_FOCE, "title", "Fraction ocean"), fnam)
275    CALL ncerr(nf90_put_att(nid, id_FSIC, "title", "Fraction glace de mer"), fnam)
276    CALL ncerr(nf90_put_att(nid, id_FTER, "title", "Fraction terre"), fnam)
277    CALL ncerr(nf90_put_att(nid, id_FLIC, "title", "Fraction land ice"), fnam)
278    CALL ncerr(nf90_put_att(nid, id_SST, "title", "Temperature superficielle de la mer"), fnam)
279    CALL ncerr(nf90_put_att(nid, id_BILS, "title", "Reference flux de chaleur au sol"), fnam)
280    CALL ncerr(nf90_put_att(nid, id_ALB, "title", "Albedo a la surface"), fnam)
281    CALL ncerr(nf90_put_att(nid, id_RUG, "title", "Rugosite"), fnam)
282
283    CALL nf95_put_att(nid, varid_longitude, "standard_name", "longitude")
284    CALL nf95_put_att(nid, varid_longitude, "units", "degrees_east")
285
286    CALL nf95_put_att(nid, varid_latitude, "standard_name", "latitude")
287    CALL nf95_put_att(nid, varid_latitude, "units", "degrees_north")
288
289    CALL ncerr(nf90_enddef(nid), fnam)
290
291    !--- Variables saving
292    CALL ncerr(nf90_put_var(nid, id_tim, [(REAL(k), k = 1, ndays)]), fnam)
293    CALL ncerr(nf90_put_var(nid, id_FOCE, pctsrf_t(:, is_oce, :), [1, 1], [klon, ndays]), fnam)
294    CALL ncerr(nf90_put_var(nid, id_FSIC, pctsrf_t(:, is_sic, :), [1, 1], [klon, ndays]), fnam)
295    CALL ncerr(nf90_put_var(nid, id_FTER, pctsrf_t(:, is_ter, :), [1, 1], [klon, ndays]), fnam)
296    CALL ncerr(nf90_put_var(nid, id_FLIC, pctsrf_t(:, is_lic, :), [1, 1], [klon, ndays]), fnam)
297    CALL ncerr(nf90_put_var(nid, id_SST, phy_sst(:, :), [1, 1], [klon, ndays]), fnam)
298    CALL ncerr(nf90_put_var(nid, id_BILS, phy_bil(:, :), [1, 1], [klon, ndays]), fnam)
299    CALL ncerr(nf90_put_var(nid, id_ALB, phy_alb(:, :), [1, 1], [klon, ndays]), fnam)
300    CALL ncerr(nf90_put_var(nid, id_RUG, phy_rug(:, :), [1, 1], [klon, ndays]), fnam)
301    CALL nf95_put_var(nid, varid_longitude, longitude_deg)
302    CALL nf95_put_var(nid, varid_latitude, latitude_deg)
303
304    CALL ncerr(nf90_close(nid), fnam)
305
306    CALL msg(0, ""); CALL msg(0, ' *** Ecriture du fichier limit : fin ***')
307
308    DEALLOCATE(pctsrf_t, phy_sst, phy_bil, phy_alb, phy_rug)
309
310
311    !===============================================================================
312
313  CONTAINS
314
315    !===============================================================================
316
317
318    !-------------------------------------------------------------------------------
319
320    SUBROUTINE get_2Dfield(fnam, varname, mode, ndays, champo, flag, mask)
321
322      !-----------------------------------------------------------------------------
323      ! Comments:
324      !   There are two assumptions concerning the NetCDF files, that are satisfied
325      !   with files that are conforming NC convention:
326      !     1) The last dimension of the variables used is the time record.
327      !     2) Dimensional variables have the same names as corresponding dimensions.
328      !-----------------------------------------------------------------------------
329      USE netcdf, ONLY: nf90_open, nf90_inq_varid, nf90_inquire_variable, &
330              nf90_close, nf90_inq_dimid, nf90_inquire_dimension, nf90_get_var, &
331              nf90_get_att
332      USE lmdz_libmath_pch, ONLY: pchsp_95, pchfe_95
333      USE lmdz_arth, ONLY: arth
334      USE indice_sol_mod
335      USE lmdz_abort_physic, ONLY: abort_physic
336      USE lmdz_libmath, ONLY: minmax
337      USE lmdz_comgeom2
338
339USE lmdz_dimensions, ONLY: iim, jjm, llm, ndm
340  USE lmdz_paramet
341      IMPLICIT NONE
342
343
344      !-----------------------------------------------------------------------------
345      ! Arguments:
346      CHARACTER(LEN = *), INTENT(IN) :: fnam     ! NetCDF file name
347      CHARACTER(LEN = *), INTENT(IN) :: varname  ! NetCDF variable name
348      CHARACTER(LEN = *), INTENT(IN) :: mode     ! RUG, SIC, SST or ALB
349      INTEGER, INTENT(IN) :: ndays    ! current year number of days
350      REAL, POINTER, DIMENSION(:, :) :: champo  ! output field = f(t)
351      LOGICAL, OPTIONAL, INTENT(IN) :: flag     ! extrapol. (SST) old ice (SIC)
352      REAL, OPTIONAL, DIMENSION(iim, jjp1), INTENT(IN) :: mask
353      !------------------------------------------------------------------------------
354      ! Local variables:
355      !--- NetCDF
356      INTEGER :: ncid, varid        ! NetCDF identifiers
357      CHARACTER(LEN = ns) :: dnam               ! dimension name
358      !--- dimensions
359      INTEGER :: dids(4)            ! NetCDF dimensions identifiers
360      REAL, ALLOCATABLE :: dlon_ini(:)        ! initial longitudes vector
361      REAL, ALLOCATABLE :: dlat_ini(:)        ! initial latitudes  vector
362      REAL, POINTER :: dlon(:), dlat(:)   ! reordered lon/lat  vectors
363      !--- fields
364      INTEGER :: imdep, jmdep, lmdep          ! dimensions of 'champ'
365      REAL, ALLOCATABLE :: champ(:, :)         ! wanted field on initial grid
366      REAL, ALLOCATABLE :: yder(:), timeyear(:)
367      REAL :: champint(iim, jjp1) ! interpolated field
368      REAL, ALLOCATABLE :: champtime(:, :, :)
369      REAL, ALLOCATABLE :: champan(:, :, :)
370      !--- input files
371      CHARACTER(LEN = ns) :: fnam_m, fnam_p     ! previous/next files names
372      CHARACTER(LEN = ns) :: cal_in             ! calendar
373      CHARACTER(LEN = ns) :: units              ! attribute "units" in sic/sst file
374      INTEGER :: ndays_in           ! number of days
375      REAL :: value              ! mean/max value near equator
376      !--- misc
377      INTEGER :: i, j, k, l         ! loop counters
378      REAL, ALLOCATABLE :: work(:, :)          ! used for extrapolation
379      CHARACTER(LEN = ns) :: title, mess        ! for messages
380      LOGICAL :: is_bcs             ! flag for BCS data
381      LOGICAL :: extrp              ! flag for extrapolation
382      LOGICAL :: ll
383      REAL :: chmin, chmax, timeday, al
384      INTEGER ierr, idx
385      INTEGER n_extrap ! number of extrapolated points
386      LOGICAL skip
387
388      !------------------------------------------------------------------------------
389      !---Variables depending on keyword 'mode' -------------------------------------
390      NULLIFY(champo)
391
392      SELECT CASE(mode)
393      CASE('RUG'); title = 'Rugosite'
394      CASE('SIC'); title = 'Sea-ice'
395      CASE('SST'); title = 'SST'
396      CASE('ALB'); title = 'Albedo'
397      END SELECT
398      extrp = .FALSE.; IF(PRESENT(flag).AND.mode=='SST') extrp = flag
399      is_bcs = (mode=='SIC'.AND.ix_sic==1).OR.(mode=='SST'.AND.ix_sst==1)
400      idx = INDEX(fnam, '.nc') - 1
401
402      !--- GETTING SOME DIMENSIONAL VARIABLES FROM FILE -----------------------------
403      CALL msg(5, ' Now reading file : ' // TRIM(fnam))
404      CALL ncerr(nf90_open(fnam, nf90_nowrite, ncid), fnam)
405      CALL ncerr(nf90_inq_varid(ncid, trim(varname), varid), fnam)
406      CALL ncerr(nf90_inquire_variable(ncid, varid, dimids = dids), fnam)
407
408      !--- Longitude
409      CALL ncerr(nf90_inquire_dimension(ncid, dids(1), name = dnam, len = imdep), fnam)
410      ALLOCATE(dlon_ini(imdep), dlon(imdep))
411      CALL ncerr(nf90_inq_varid(ncid, dnam, varid), fnam)
412      CALL ncerr(nf90_get_var(ncid, varid, dlon_ini), fnam)
413      CALL msg(5, 'variable ' // TRIM(dnam) // ' dimension ', imdep)
414
415      !--- Latitude
416      CALL ncerr(nf90_inquire_dimension(ncid, dids(2), name = dnam, len = jmdep), fnam)
417      ALLOCATE(dlat_ini(jmdep), dlat(jmdep))
418      CALL ncerr(nf90_inq_varid(ncid, dnam, varid), fnam)
419      CALL ncerr(nf90_get_var(ncid, varid, dlat_ini), fnam)
420      CALL msg(5, 'variable ' // TRIM(dnam) // ' dimension ', jmdep)
421
422      !--- Time (variable is not needed - it is rebuilt - but calendar is)
423      CALL ncerr(nf90_inquire_dimension(ncid, dids(3), name = dnam, len = lmdep), fnam)
424      ALLOCATE(timeyear(lmdep + 2))
425      CALL ncerr(nf90_inq_varid(ncid, dnam, varid), fnam)
426      cal_in = ' '
427      IF(nf90_get_att(ncid, varid, 'calendar', cal_in)/=nf90_noerr) THEN
428        SELECT CASE(mode)
429        CASE('RUG', 'ALB'); cal_in = '360_day'
430        CASE('SIC', 'SST'); cal_in = 'gregorian'
431        END SELECT
432        CALL msg(0, 'WARNING: missing "calendar" attribute for "time" in '&
433                // TRIM(fnam) // '. Choosing default value.')
434      END IF
435      CALL strclean(cal_in)                     !--- REMOVE (WEIRD) NULL CHARACTERS
436      CALL msg(0, 'var, calendar, dim: ' // TRIM(dnam) // ' ' // TRIM(cal_in), lmdep)
437
438      !--- Determining input file number of days, depending on calendar
439      ndays_in = year_len(anneeref, cal_in)
440
441      !--- Rebuilding input time vector (field from input file might be unreliable)
442      IF(lmdep==12) THEN
443        timeyear = mid_month(anneeref, cal_in)
444        CALL msg(0, 'Monthly input file(s) for ' // TRIM(title) // '.')
445      ELSE IF(lmdep==ndays_in) THEN
446        timeyear = [(REAL(k) - 0.5, k = 0, ndays_in + 1)]
447        CALL msg(0, 'Daily input file (no time interpolation).')
448      ELSE
449        WRITE(mess, '(a,i3,a,i3,a)')'Mismatching input file: found', lmdep, &
450                ' records, 12/', ndays_in, ' (monthly/daily needed).'
451        CALL abort_physic('mid_month', TRIM(mess), 1)
452      END IF
453
454      !--- GETTING THE FIELD AND INTERPOLATING IT ----------------------------------
455      ALLOCATE(champ(imdep, jmdep), champtime(iim, jjp1, lmdep + 2))
456      IF(extrp) ALLOCATE(work(imdep, jmdep))
457      CALL msg(5, '')
458      CALL msg(5, 'READ AND INTERPOLATE HORIZONTALLY ', lmdep, ' FIELDS.')
459      CALL ncerr(nf90_inq_varid(ncid, varname, varid), fnam)
460      DO l = 1, lmdep
461        CALL ncerr(nf90_get_var(ncid, varid, champ, [1, 1, l], [imdep, jmdep, 1]), fnam)
462        CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.)
463
464        !--- FOR SIC/SST FIELDS ONLY
465        IF(l==1.AND.is_in(mode, ['SIC', 'SST'])) THEN
466
467          !--- DETERMINE THE UNIT: READ FROM FILE OR ASSUMED USING FIELD VALUES
468          ierr = nf90_get_att(ncid, varid, 'units', units)
469          IF(ierr==nf90_noerr) THEN !--- ATTRIBUTE "units" FOUND IN THE FILE
470            CALL strclean(units)
471            IF(mode=='SIC'.AND.is_in(units, Perc)) units = "%"
472            IF(mode=='SIC'.AND.is_in(units, Frac)) units = "1"
473            IF(mode=='SST'.AND.is_in(units, DegC)) units = "C"
474            IF(mode=='SST'.AND.is_in(units, DegK)) units = "K"
475          ELSE                      !--- CHECK THE FIELD VALUES
476            IF(mode=='SIC') value = MAXVAL(champ(:, :))
477            IF(mode=='SST') value = SUM(champ(:, jmdep / 2), DIM = 1) / REAL(imdep)
478            IF(mode=='SIC') THEN; units = "1"; IF(value>= 10.) units = "%";
479            END IF
480            IF(mode=='SST') THEN; units = "C"; IF(value>=100.) units = "K";
481            END IF
482          END IF
483          CALL msg(0, 'INPUT FILE ' // TRIM(title) // ' UNIT IS: "' // TRIM(units) // '".')
484          IF(ierr/=nf90_noerr) CALL msg(0, 'WARNING ! UNIT TO BE CHECKED ! '      &
485                  // 'No "units" attribute, so only based on the fields values.')
486
487          !--- CHECK VALUES ARE IN THE EXPECTED RANGE
488          SELECT CASE(units)
489          CASE('%'); ll = ANY(champ>100.0 + EPSFRA); str = 'percentages > 100.'
490          CASE('1'); ll = ANY(champ>  1.0 + EPSFRA); str = 'fractions > 1.'
491          CASE('C'); ll = ANY(champ<-100.).OR.ANY(champ> 60.); str = '<-100 or >60 DegC'
492          CASE('K'); ll = ANY(champ< 180.).OR.ANY(champ>330.); str = '<180 or >330 DegK'
493          CASE DEFAULT; CALL abort_physic(mode, 'Unrecognized ' // TRIM(title)   &
494                  // ' unit: ' // TRIM(units), 1)
495          END SELECT
496
497          !--- DROPPED FOR BCS DATA (FRACTIONS CAN BE HIGHER THAN 1)
498          IF(ll.AND.ix_sic/=1.AND.mode=='SIC') &
499                  CALL abort_physic(mode, 'unrealistic ' // TRIM(mode) // ' found: ' // TRIM(str), 1)
500
501        END IF
502
503        IF(extrp) CALL extrapol(champ, imdep, jmdep, 999999., .TRUE., .TRUE., 2, work)
504        IF(l==1) THEN
505          CALL msg(5, "--------------------------------------------------------")
506          CALL msg(5, "$$$ Barycentric interpolation for " // TRIM(title) // " $$$")
507          CALL msg(5, "--------------------------------------------------------")
508        END IF
509        IF(mode=='RUG') champ = LOG(champ)
510        CALL inter_barxy(dlon, dlat(:jmdep - 1), champ, rlonu(:iim), rlatv, champint)
511        IF(mode=='RUG') THEN
512          champint = EXP(champint)
513          WHERE(NINT(mask)/=1) champint = 0.001
514        END IF
515        champtime(:, :, l + 1) = champint
516      END DO
517      CALL ncerr(nf90_close(ncid), fnam)
518
519      !--- FIRST RECORD: LAST ONE OF PREVIOUS YEAR (CURRENT YEAR IF UNAVAILABLE)
520      fnam_m = fnam(1:idx) // '_m.nc'
521      IF(nf90_open(fnam_m, nf90_nowrite, ncid)==nf90_noerr) THEN
522        CALL msg(0, 'Reading previous year file ("' // TRIM(fnam_m) // '") last record for ' // TRIM(title))
523        CALL ncerr(nf90_inq_varid(ncid, varname, varid), fnam_m)
524        CALL ncerr(nf90_inquire_variable(ncid, varid, dimids = dids), fnam_m)
525        CALL ncerr(nf90_inquire_dimension(ncid, dids(3), len = l), fnam_m)
526        CALL ncerr(nf90_get_var(ncid, varid, champ, [1, 1, l], [imdep, jmdep, 1]), fnam_m)
527        CALL ncerr(nf90_close(ncid), fnam_m)
528        CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.)
529        IF(extrp) CALL extrapol(champ, imdep, jmdep, 999999., .TRUE., .TRUE., 2, work)
530        IF(mode=='RUG') champ = LOG(champ)
531        CALL inter_barxy(dlon, dlat(:jmdep - 1), champ, rlonu(:iim), rlatv, champint)
532        IF(mode=='RUG') THEN
533          champint = EXP(champint)
534          WHERE(NINT(mask)/=1) champint = 0.001
535        END IF
536        champtime(:, :, 1) = champint
537      ELSE
538        CALL msg(0, 'Using current year file ("' // TRIM(fnam) // '") last record for ' // TRIM(title))
539        champtime(:, :, 1) = champtime(:, :, lmdep + 1)
540      END IF
541
542      !--- LAST RECORD: FIRST ONE OF NEXT YEAR (CURRENT YEAR IF UNAVAILABLE)
543      fnam_p = fnam(1:idx) // '_p.nc'
544      IF(nf90_open(fnam_p, nf90_nowrite, ncid)==nf90_noerr) THEN
545        CALL msg(0, 'Reading next year file ("' // TRIM(fnam_p) // '") first record for ' // TRIM(title))
546        CALL ncerr(nf90_inq_varid(ncid, varname, varid), fnam_p)
547        CALL ncerr(nf90_get_var(ncid, varid, champ, [1, 1, 1], [imdep, jmdep, 1]), fnam_p)
548        CALL ncerr(nf90_close(ncid), fnam_p)
549        CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.)
550        IF(extrp) CALL extrapol(champ, imdep, jmdep, 999999., .TRUE., .TRUE., 2, work)
551        IF(mode=='RUG') champ = LOG(champ)
552        CALL inter_barxy(dlon, dlat(:jmdep - 1), champ, rlonu(:iim), rlatv, champint)
553        IF(mode=='RUG') THEN
554          champint = EXP(champint)
555          WHERE(NINT(mask)/=1) champint = 0.001
556        END IF
557        champtime(:, :, lmdep + 2) = champint
558      ELSE
559        CALL msg(0, 'Using current year file ("' // TRIM(fnam) // '") first record for ' // TRIM(title))
560        champtime(:, :, lmdep + 2) = champtime(:, :, 2)
561      END IF
562      DEALLOCATE(dlon_ini, dlat_ini, dlon, dlat, champ)
563      IF(extrp) DEALLOCATE(work)
564
565      !--- TIME INTERPOLATION ------------------------------------------------------
566      IF(prt_level>0) THEN
567        IF(ndays/=ndays_in) THEN
568          WRITE(lunout, *)'DIFFERENT YEAR LENGTHS:'
569          WRITE(lunout, *)' In the  input file: ', ndays_in
570          WRITE(lunout, *)' In the output file: ', ndays
571        END IF
572        IF(lmdep==ndays_in) THEN
573          WRITE(lunout, *)'NO TIME INTERPOLATION.'
574          WRITE(lunout, *)' Daily input file.'
575        ELSE
576          IF(is_bcs) WRITE(lunout, *)'LINEAR TIME INTERPOLATION.'
577          IF(.NOT.is_bcs) WRITE(lunout, *)'SPLINES TIME INTERPOLATION.'
578          WRITE(lunout, *)' Input time vector: ', timeyear
579          WRITE(lunout, *)' Output time vector: from 0.5 to ', ndays - 0.5
580        END IF
581      END IF
582      ALLOCATE(champan(iip1, jjp1, ndays))
583
584      IF(lmdep==ndays_in) THEN  !--- DAILY DATA: NO     TIME INTERPOLATION
585        DO l = 1, lmdep
586          champan(1:iim, :, l) = champtime(:, :, l + 1)
587        END DO
588      ELSE IF(is_bcs) THEN      !--- BCS   DATA: LINEAR TIME INTERPOLATION
589        l = 1
590        DO k = 1, ndays
591          timeday = (REAL(k) - 0.5) * REAL(ndays_in) / ndays
592          IF(timeyear(l + 1)<timeday) l = l + 1
593          al = (timeday - timeyear(l)) / (timeyear(l + 1) - timeyear(l))
594          champan(1:iim, :, k) = champtime(1:iim, :, l) + al * (champtime(1:iim, :, l + 1) - champtime(1:iim, :, l))
595        END DO
596      ELSE                      !--- AVE   DATA: SPLINE TIME INTERPOLATION
597        skip = .FALSE.
598        n_extrap = 0
599        ALLOCATE(yder(lmdep + 2))
600        DO j = 1, jjp1
601          DO i = 1, iim
602            yder = pchsp_95(timeyear, champtime(i, j, :), ibeg = 2, iend = 2, &
603                    vc_beg = 0., vc_end = 0.)
604            CALL pchfe_95(timeyear, champtime(i, j, :), yder, skip, &
605                    arth(0.5, real(ndays_in) / ndays, ndays), champan(i, j, :), ierr)
606            IF (ierr < 0) CALL abort_physic("get_2Dfield", "", 1)
607            n_extrap = n_extrap + ierr
608          END DO
609        END DO
610        IF(n_extrap /= 0) WRITE(lunout, *) "get_2Dfield pchfe_95: n_extrap = ", n_extrap
611        DEALLOCATE(yder)
612      END IF
613      champan(iip1, :, :) = champan(1, :, :)
614      DEALLOCATE(champtime, timeyear)
615
616      !--- Checking the result
617      DO j = 1, jjp1
618        CALL minmax(iip1, champan(1, j, 10), chmin, chmax)
619        IF (prt_level>5) WRITE(lunout, *)' ', TRIM(title), ' at time 10 ', chmin, chmax, j
620      END DO
621
622      !--- SPECIAL FILTER FOR SST: SST>271.38 --------------------------------------
623      IF(mode=='SST') THEN
624        SELECT CASE(units)
625        CASE("K"); CALL msg(0, 'SST field is already in kelvins.')
626        CASE("C"); CALL msg(0, 'SST field converted from celcius degrees to kelvins.')
627        champan(:, :, :) = champan(:, :, :) + 273.15
628        END SELECT
629        CALL msg(0, 'Filtering SST: Sea Surface Temperature >= 271.38')
630        WHERE(champan<271.38) champan = 271.38
631      END IF
632
633      !--- SPECIAL FILTER FOR SIC: 0.0<SIC<1.0 -------------------------------------
634      IF(mode=='SIC') THEN
635        SELECT CASE(units)
636        CASE("1"); CALL msg(0, 'SIC field already in fraction of 1')
637        CASE("%"); CALL msg(0, 'SIC field converted from percentage to fraction of 1.')
638        champan(:, :, :) = champan(:, :, :) / 100.
639        END SELECT
640        CALL msg(0, 'Filtering SIC: 0.0 <= Sea-ice <=1.0')
641        WHERE(champan>1.0) champan = 1.0
642        WHERE(champan<0.0) champan = 0.0
643      END IF
644
645      !--- DYNAMICAL TO PHYSICAL GRID ----------------------------------------------
646      ALLOCATE(champo(klon, ndays))
647      DO k = 1, ndays
648        CALL gr_dyn_fi(1, iip1, jjp1, klon, champan(1, 1, k), champo(:, k))
649      END DO
650      DEALLOCATE(champan)
651
652    END SUBROUTINE get_2Dfield
653
654    !-------------------------------------------------------------------------------
655
656
657    !-------------------------------------------------------------------------------
658
659    SUBROUTINE start_init_orog0(lon_in, lat_in, phis, masque)
660      USE grid_noro_m, ONLY: grid_noro0
661
662
663      IMPLICIT NONE
664      !===============================================================================
665      ! Purpose:  Compute "phis" just like it would be in start_init_orog.
666      !===============================================================================
667      ! Arguments:
668      REAL, INTENT(IN) :: lon_in(:), lat_in(:)   ! dim (iml) (jml)
669      REAL, INTENT(INOUT) :: phis(:, :), masque(:, :) ! dim (iml,jml)
670      !-------------------------------------------------------------------------------
671      ! Local variables:
672      CHARACTER(LEN = ns) :: modname = "start_init_orog0"
673      INTEGER :: fid, llm_tmp, ttm_tmp, iml, jml, iml_rel, jml_rel, itau(1)
674      REAL :: lev(1), date, dt, deg2rad
675      REAL, ALLOCATABLE :: lon_rad(:), lon_ini(:), lon_rel(:, :), relief_hi(:, :)
676      REAL, ALLOCATABLE :: lat_rad(:), lat_ini(:), lat_rel(:, :)
677      !-------------------------------------------------------------------------------
678      iml = assert_eq(SIZE(lon_in), SIZE(phis, 1), SIZE(masque, 1), TRIM(modname) // " iml")
679      jml = assert_eq(SIZE(lat_in), SIZE(phis, 2), SIZE(masque, 2), TRIM(modname) // " jml")
680      IF(iml/=iip1) CALL abort_gcm(TRIM(modname), 'iml/=iip1', 1)
681      IF(jml/=jjp1) CALL abort_gcm(TRIM(modname), 'jml/=jjp1', 1)
682      pi = 2.0 * ASIN(1.0); deg2rad = pi / 180.0
683      IF(ANY(phis/=-99999.)) RETURN                  !--- phis ALREADY KNOWN
684
685      !--- HIGH RESOLUTION OROGRAPHY
686      CALL flininfo(frelf, iml_rel, jml_rel, llm_tmp, ttm_tmp, fid)
687
688      ALLOCATE(lat_rel(iml_rel, jml_rel), lon_rel(iml_rel, jml_rel))
689      CALL flinopen(frelf, .FALSE., iml_rel, jml_rel, llm_tmp, lon_rel, lat_rel, &
690              lev, ttm_tmp, itau, date, dt, fid)
691      ALLOCATE(relief_hi(iml_rel, jml_rel))
692      CALL flinget(fid, vrel, iml_rel, jml_rel, llm_tmp, ttm_tmp, 1, 1, relief_hi)
693      CALL flinclo(fid)
694
695      !--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS
696      ALLOCATE(lon_ini(iml_rel), lat_ini(jml_rel))
697      lon_ini(:) = lon_rel(:, 1); IF(MAXVAL(lon_rel)>pi) lon_ini = lon_ini * deg2rad
698      lat_ini(:) = lat_rel(1, :); IF(MAXVAL(lat_rel)>pi) lat_ini = lat_ini * deg2rad
699
700      !--- FIELDS ARE PROCESSED TO BE ON STANDARD ANGULAR DOMAINS
701      ALLOCATE(lon_rad(iml_rel), lat_rad(jml_rel))
702      CALL conf_dat2d(vrel, lon_ini, lat_ini, lon_rad, lat_rad, relief_hi, .FALSE.)
703      DEALLOCATE(lon_ini, lat_ini)
704
705      !--- COMPUTING SURFACE GEOPOTENTIAL USING ROUTINE grid_noro0
706      WRITE(lunout, *)
707      WRITE(lunout, *)'*** Compute surface geopotential ***'
708
709      !--- CALL OROGRAPHY MODULE (REDUCED VERSION) TO COMPUTE FIELDS
710      CALL grid_noro0(lon_rad, lat_rad, relief_hi, lon_in, lat_in, phis, masque)
711      phis = phis * 9.81
712      phis(iml, :) = phis(1, :)
713      DEALLOCATE(relief_hi, lon_rad, lat_rad)
714
715    END SUBROUTINE start_init_orog0
716
717    !-------------------------------------------------------------------------------
718
719
720    !-------------------------------------------------------------------------------
721
722    SUBROUTINE msg(lev, str1, i, str2)
723
724      !-------------------------------------------------------------------------------
725      ! Arguments:
726      INTEGER, INTENT(IN) :: lev
727      CHARACTER(LEN = *), INTENT(IN) :: str1
728      INTEGER, OPTIONAL, INTENT(IN) :: i
729      CHARACTER(LEN = *), OPTIONAL, INTENT(IN) :: str2
730      !-------------------------------------------------------------------------------
731      IF(prt_level>=lev) THEN
732        IF(PRESENT(str2)) THEN
733          WRITE(lunout, *) TRIM(str1), i, TRIM(str2)
734        ELSE IF(PRESENT(i)) THEN
735          WRITE(lunout, *) TRIM(str1), i
736        ELSE
737          WRITE(lunout, *) TRIM(str1)
738        END IF
739      END IF
740
741    END SUBROUTINE msg
742
743    !-------------------------------------------------------------------------------
744
745
746    !-------------------------------------------------------------------------------
747
748    SUBROUTINE ncerr(ncres, fnam)
749
750      !-------------------------------------------------------------------------------
751      ! Purpose: NetCDF errors handling.
752      !-------------------------------------------------------------------------------
753      USE netcdf, ONLY: nf90_noerr, nf90_strerror
754      USE lmdz_abort_physic, ONLY: abort_physic
755      IMPLICIT NONE
756      !-------------------------------------------------------------------------------
757      ! Arguments:
758      INTEGER, INTENT(IN) :: ncres
759      CHARACTER(LEN = *), INTENT(IN) :: fnam
760      !-------------------------------------------------------------------------------
761      IF(ncres/=nf90_noerr) THEN
762        WRITE(lunout, *)'Problem with file ' // TRIM(fnam) // ' in routine limit_netcdf.'
763        CALL abort_physic('limit_netcdf', nf90_strerror(ncres), 1)
764      END IF
765
766    END SUBROUTINE ncerr
767
768    !-------------------------------------------------------------------------------
769
770
771    !-------------------------------------------------------------------------------
772
773    SUBROUTINE strclean(s)
774
775      !-------------------------------------------------------------------------------
776      IMPLICIT NONE
777      !-------------------------------------------------------------------------------
778      ! Purpose: Remove tail null characters from the input string.
779      !-------------------------------------------------------------------------------
780      ! Parameters:
781      CHARACTER(LEN = *), INTENT(INOUT) :: s
782      !-------------------------------------------------------------------------------
783      ! Local variable:
784      INTEGER :: k
785      !-------------------------------------------------------------------------------
786      k = LEN_TRIM(s); DO WHILE(ICHAR(s(k:k))==0); s(k:k) = ' '; k = LEN_TRIM(s);
787      END DO
788
789    END SUBROUTINE strclean
790
791    !-------------------------------------------------------------------------------
792
793
794    !-------------------------------------------------------------------------------
795
796    FUNCTION is_in(s1, s2) RESULT(res)
797
798      !-------------------------------------------------------------------------------
799      IMPLICIT NONE
800      !-------------------------------------------------------------------------------
801      ! Purpose: Check wether s1 is present in the s2(:) list (case insensitive).
802      !-------------------------------------------------------------------------------
803      ! Arguments:
804      CHARACTER(LEN = *), INTENT(IN) :: s1, s2(:)
805      LOGICAL :: res
806      !-------------------------------------------------------------------------------
807      res = .FALSE.; DO k = 1, SIZE(s2); res = res.OR.strLow(s1)==strLow(s2(k));
808      END DO
809
810    END FUNCTION is_in
811
812    !-------------------------------------------------------------------------------
813
814
815    !-------------------------------------------------------------------------------
816
817    ELEMENTAL FUNCTION strLow(s) RESULT(res)
818
819      !-------------------------------------------------------------------------------
820      IMPLICIT NONE
821      !-------------------------------------------------------------------------------
822      ! Purpose: Lower case conversion.
823      !-------------------------------------------------------------------------------
824      ! Arguments:
825      CHARACTER(LEN = *), INTENT(IN) :: s
826      CHARACTER(LEN = ns) :: res
827      !-------------------------------------------------------------------------------
828      ! Local variable:
829      INTEGER :: k, ix
830      !-------------------------------------------------------------------------------
831      res = s
832      DO k = 1, LEN(s); ix = IACHAR(s(k:k))
833      IF(64<ix.AND.ix<91) res(k:k) = ACHAR(ix + 32)
834      END DO
835
836    END FUNCTION strLow
837
838    !-------------------------------------------------------------------------------
839
840  END SUBROUTINE limit_netcdf
841
842END MODULE limit
843
844!*******************************************************************************
845
Note: See TracBrowser for help on using the repository browser.