source: LMDZ6/trunk/libf/dynphy_lonlat/phylmd/limit_netcdf.F90 @ 5073

Last change on this file since 5073 was 5073, checked in by abarral, 7 weeks ago

Remove all NC_DOUBLE uses outside of lmdz_netcdf.F90 (except in obsolete/, which I hope we'll ditch soon...)
Note: make sure to check convergence at some point, it's possible that we've messed up some when replacing nf_* by nf90_* calls
(lint) replace obsolete logical operators along the way

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