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

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

Turn paramet.h into a module

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