GCC Code Coverage Report
Directory: ./ Exec Total Coverage
File: dynphy_lonlat/phylmd/limit_netcdf.F90 Lines: 0 360 0.0 %
Date: 2023-06-30 12:51:15 Branches: 0 922 0.0 %

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