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

Last change on this file since 5113 was 5113, checked in by abarral, 2 months ago

Rename modules in misc from *_mod > lmdz_*
Put cbrt.f90, ch*.f90, pch*.f90 in new lmdz_libmath_pch.f90

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