source: LMDZ6/branches/Amaury_dev/libf/dynphy_lonlat/phylmd/limit_netcdf.F90 @ 5080

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

[continued & end] replace netcdf by lmdz_netcdf.F90 wrapper
"use netcdf" is now only used in lmdz_netcdf.F90 (except ecrad and obsolete/)
<include "netcdf.inc"> is now likewise only used in lmdz_netcdf.F90.

systematically specify explicitely <USE lmdz_netcdf, ONLY:> (probably left some missing, to correct later on)

Further replacement of nf_put_* by nf90_put_* (same for _get_)

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