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

Last change on this file since 5139 was 5136, checked in by abarral, 8 weeks ago

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