source: LMDZ5/branches/testing/libf/dynphy_lonlat/phylmd/limit_netcdf.F90 @ 2542

Last change on this file since 2542 was 2542, checked in by Laurent Fairhead, 8 years ago

Merged trunk changes r2487:2541 into testing branch

  • 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.4 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   ioconf_calendar, ioget_calendar, lock_calendar, ioget_mon_len, ioget_year_len
21  USE assert_eq_m,        ONLY: assert_eq
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  CHARACTER(LEN=20), PARAMETER :: &
30  fsst(4)=['amipbc_sst_1x1.nc   ','cpl_atm_sst.nc      ','histmth_sst.nc      '&
31          ,'sstk.nc             ']
32  CHARACTER(LEN=20), PARAMETER :: &
33  fsic(4)=['amipbc_sic_1x1.nc   ','cpl_atm_sic.nc      ','histmth_sic.nc      '&
34          ,'ci.nc               ']
35  CHARACTER(LEN=10), PARAMETER :: &
36  vsst(4)=['tosbcs    ','SISUTESW  ','tsol_oce  ','sstk      '], &
37  vsic(4)=['sicbcs    ','SIICECOV  ','pourc_sic ','ci        ']
38  CHARACTER(LEN=10), PARAMETER :: &
39  frugo='Rugos.nc  ', falbe='Albedo.nc ', frelf='Relief.nc ',    &
40   vrug='RUGOS     ',  valb='ALBEDO    ',  vrel='RELIEF    '
41
42CONTAINS
43
44!-------------------------------------------------------------------------------
45!
46SUBROUTINE limit_netcdf(masque, phis, extrap)
47!
48!-------------------------------------------------------------------------------
49! Author : L. Fairhead, 27/01/94
50!-------------------------------------------------------------------------------
51! Purpose: Boundary conditions files building for new model using climatologies.
52!          Both grids have to be regular.
53!-------------------------------------------------------------------------------
54! Note: This routine is designed to work for Earth
55!-------------------------------------------------------------------------------
56! Modification history:
57!  * 23/03/1994: Z. X. Li
58!  *    09/1999: L. Fairhead (netcdf reading in LMDZ.3.3)
59!  *    07/2001: P. Le Van
60!  *    11/2009: L. Guez     (ozone day & night climatos, see etat0_netcdf.F90)
61!  *    12/2009: D. Cugnet   (f77->f90, calendars, files from coupled runs)
62!  *    04/2016: D. Cugnet   (12/14 recs SST/SIC files: cyclic/interannual runs)
63!-------------------------------------------------------------------------------
64#ifndef CPP_1D
65  USE indice_sol_mod
66  USE netcdf,             ONLY: NF90_OPEN,    NF90_CREATE,  NF90_CLOSE,        &
67                  NF90_DEF_DIM, NF90_DEF_VAR, NF90_PUT_VAR, NF90_PUT_ATT,      &
68                  NF90_NOERR,   NF90_NOWRITE, NF90_DOUBLE,  NF90_GLOBAL,       &
69                  NF90_CLOBBER, NF90_ENDDEF,  NF90_UNLIMITED, NF90_FLOAT
70  USE inter_barxy_m,      ONLY: inter_barxy
71  USE netcdf95,           ONLY: nf95_def_var, nf95_put_att, nf95_put_var
72  IMPLICIT NONE
73!-------------------------------------------------------------------------------
74! Arguments:
75  include "iniprint.h"
76  include "dimensions.h"
77  include "paramet.h"
78  REAL, DIMENSION(iip1,jjp1), INTENT(INOUT) :: masque ! land mask
79  REAL, DIMENSION(iip1,jjp1), INTENT(INOUT) :: phis   ! ground geopotential
80  LOGICAL,                    INTENT(IN)    :: extrap ! SST extrapolation flag
81!-------------------------------------------------------------------------------
82! Local variables:
83  include "logic.h"
84  include "comgeom2.h"
85  include "comconst.h"
86
87!--- INPUT NETCDF FILES NAMES --------------------------------------------------
88  CHARACTER(LEN=20) :: icefile, sstfile, dumstr, fnam
89  CHARACTER(LEN=10) :: varname
90
91!--- OUTPUT VARIABLES FOR NETCDF FILE ------------------------------------------
92  REAL               :: fi_ice(klon), verif(klon)
93  REAL, POINTER      :: phy_rug(:,:)=>NULL(), phy_ice(:,:)=>NULL()
94  REAL, POINTER      :: phy_sst(:,:)=>NULL(), phy_alb(:,:)=>NULL()
95  REAL, ALLOCATABLE  :: phy_bil(:,:), pctsrf_t(:,:,:)
96  INTEGER            :: nbad
97
98!--- VARIABLES FOR OUTPUT FILE WRITING -----------------------------------------
99  INTEGER :: ierr, nid, ndim, ntim, k, dims(2), ix_sic, ix_sst
100  INTEGER :: id_tim,  id_SST,  id_BILS, id_RUG, id_ALB
101  INTEGER :: id_FOCE, id_FSIC, id_FTER, id_FLIC, varid_longitude, varid_latitude
102  INTEGER :: NF90_FORMAT
103  INTEGER :: ndays                   !--- Depending on the output calendar
104
105!--- INITIALIZATIONS -----------------------------------------------------------
106#ifdef NC_DOUBLE
107  NF90_FORMAT=NF90_DOUBLE
108#else
109  NF90_FORMAT=NF90_FLOAT
110#endif
111  CALL inigeom
112
113!--- MASK, GROUND GEOPOT. & SUBSURFACES COMPUTATION (IN CASE ok_etat0==.FALSE.)
114   IF(ALL(masque==-99999.)) THEN
115    CALL start_init_orog0(rlonv,rlatu,phis,masque)
116    CALL gr_dyn_fi(1,iip1,jjp1,klon,masque,zmasq)          !--- To physical grid
117    ALLOCATE(pctsrf(klon,nbsrf))
118    CALL start_init_subsurf(.FALSE.)
119  END IF
120
121!--- Beware: anneeref (from gcm.def) is used to determine output time sampling
122  ndays=ioget_year_len(anneeref)
123
124!--- RUGOSITY TREATMENT --------------------------------------------------------
125  CALL msg(1,'Traitement de la rugosite')
126  CALL get_2Dfield(frugo,vrug,'RUG',ndays,phy_rug,mask=masque(1:iim,:))
127
128!--- OCEAN TREATMENT -----------------------------------------------------------
129  CALL msg(1,'Traitement de la glace oceanique')
130
131! Input SIC file selection
132! Open file only to test if available
133  DO ix_sic=1,SIZE(fsic)
134     IF ( NF90_OPEN(TRIM(fsic(ix_sic)),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
135        icefile=fsic(ix_sic); varname=vsic(ix_sic); EXIT
136     END IF
137  END DO
138  IF(ix_sic==SIZE(fsic)+1) THEN
139     WRITE(lunout,*) 'ERROR! No sea-ice input file was found.'
140     WRITE(lunout,*) 'One of following files must be available : '
141     DO k=1,SIZE(fsic); WRITE(lunout,*) TRIM(fsic(k)); END DO
142     CALL abort_physic('limit_netcdf','No sea-ice file was found',1)
143  END IF
144  CALL ncerr(NF90_CLOSE(nid),icefile)
145  CALL msg(-1,'Fichier choisi pour la glace de mer:'//TRIM(icefile))
146
147  CALL get_2Dfield(icefile,varname, 'SIC',ndays,phy_ice)
148
149  ALLOCATE(pctsrf_t(klon,nbsrf,ndays))
150  DO k=1,ndays
151     fi_ice=phy_ice(:,k)
152     WHERE(fi_ice>=1.0  ) fi_ice=1.0
153     WHERE(fi_ice<EPSFRA) fi_ice=0.0
154     pctsrf_t(:,is_ter,k)=pctsrf(:,is_ter)       ! land soil
155     pctsrf_t(:,is_lic,k)=pctsrf(:,is_lic)       ! land ice
156     SELECT CASE(ix_sic)
157        CASE(2)                                   ! SIC=pICE*(1-LIC-TER) (CPL)
158        pctsrf_t(:,is_sic,k)=fi_ice(:)*(1.-pctsrf(:,is_lic)-pctsrf(:,is_ter))
159        CASE(3)                                   ! SIC=pICE            (HIST)
160        pctsrf_t(:,is_sic,k)=fi_ice(:)
161        CASE DEFAULT                              ! SIC=pICE-LIC   (AMIP,ERAI)
162        pctsrf_t(:,is_sic,k)=fi_ice-pctsrf_t(:,is_lic,k)
163     END SELECT
164     WHERE(pctsrf_t(:,is_sic,k)<=0) pctsrf_t(:,is_sic,k)=0.
165     WHERE(1.0-zmasq<EPSFRA)
166        pctsrf_t(:,is_sic,k)=0.0
167        pctsrf_t(:,is_oce,k)=0.0
168     ELSEWHERE
169        WHERE(pctsrf_t(:,is_sic,k)>=1.0-zmasq)
170           pctsrf_t(:,is_sic,k)=1.0-zmasq
171           pctsrf_t(:,is_oce,k)=0.0
172        ELSEWHERE
173           pctsrf_t(:,is_oce,k)=1.0-zmasq-pctsrf_t(:,is_sic,k)
174           WHERE(pctsrf_t(:,is_oce,k)<EPSFRA)
175              pctsrf_t(:,is_oce,k)=0.0
176              pctsrf_t(:,is_sic,k)=1.0-zmasq
177           END WHERE
178        END WHERE
179     END WHERE
180     nbad=COUNT(pctsrf_t(:,is_oce,k)<0.0)
181     IF(nbad>0) WRITE(lunout,*) 'pb sous maille pour nb points = ',nbad
182     nbad=COUNT(ABS(SUM(pctsrf_t(:,:,k),DIM=2)-1.0)>EPSFRA)
183     IF(nbad>0) WRITE(lunout,*) 'pb sous surface pour nb points = ',nbad
184  END DO
185  DEALLOCATE(phy_ice)
186
187!--- SST TREATMENT -------------------------------------------------------------
188  CALL msg(1,'Traitement de la sst')
189
190! Input SST file selection
191! Open file only to test if available
192  DO ix_sst=1,SIZE(fsst)
193     IF ( NF90_OPEN(TRIM(fsst(ix_sst)),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
194       sstfile=fsst(ix_sst); varname=vsst(ix_sst); EXIT
195     END IF
196  END DO
197  IF(ix_sst==SIZE(fsst)+1) THEN
198     WRITE(lunout,*) 'ERROR! No sst input file was found.'
199     WRITE(lunout,*) 'One of following files must be available : '
200     DO k=1,SIZE(fsst); WRITE(lunout,*) TRIM(fsst(k)); END DO
201     CALL abort_physic('limit_netcdf','No sst file was found',1)
202  END IF
203  CALL ncerr(NF90_CLOSE(nid),sstfile)
204  CALL msg(-1,'Fichier choisi pour la temperature de mer: '//TRIM(sstfile))
205
206  CALL get_2Dfield(sstfile,varname,'SST',ndays,phy_sst,flag=extrap)
207
208!--- ALBEDO TREATMENT ----------------------------------------------------------
209  CALL msg(1,'Traitement de l albedo')
210  CALL get_2Dfield(falbe,valb,'ALB',ndays,phy_alb)
211
212!--- REFERENCE GROUND HEAT FLUX TREATMENT --------------------------------------
213  ALLOCATE(phy_bil(klon,ndays)); phy_bil=0.0
214
215!--- OUTPUT FILE WRITING -------------------------------------------------------
216  CALL msg(5,'Ecriture du fichier limit : debut')
217  fnam="limit.nc"
218
219  !--- File creation
220  CALL ncerr(NF90_CREATE(fnam,NF90_CLOBBER,nid),fnam)
221  CALL ncerr(NF90_PUT_ATT(nid,NF90_GLOBAL,"title","Fichier conditions aux limites"),fnam)
222
223  !--- Dimensions creation
224  CALL ncerr(NF90_DEF_DIM(nid,"points_physiques",klon,ndim),fnam)
225  CALL ncerr(NF90_DEF_DIM(nid,"time",NF90_UNLIMITED,ntim),fnam)
226
227  dims=[ndim,ntim]
228
229  !--- Variables creation
230  CALL ncerr(NF90_DEF_VAR(nid,"TEMPS",NF90_FORMAT,[ntim],id_tim),fnam)
231  CALL ncerr(NF90_DEF_VAR(nid,"FOCE", NF90_FORMAT,dims,id_FOCE),fnam)
232  CALL ncerr(NF90_DEF_VAR(nid,"FSIC", NF90_FORMAT,dims,id_FSIC),fnam)
233  CALL ncerr(NF90_DEF_VAR(nid,"FTER", NF90_FORMAT,dims,id_FTER),fnam)
234  CALL ncerr(NF90_DEF_VAR(nid,"FLIC", NF90_FORMAT,dims,id_FLIC),fnam)
235  CALL ncerr(NF90_DEF_VAR(nid,"SST",  NF90_FORMAT,dims,id_SST),fnam)
236  CALL ncerr(NF90_DEF_VAR(nid,"BILS", NF90_FORMAT,dims,id_BILS),fnam)
237  CALL ncerr(NF90_DEF_VAR(nid,"ALB",  NF90_FORMAT,dims,id_ALB),fnam)
238  CALL ncerr(NF90_DEF_VAR(nid,"RUG",  NF90_FORMAT,dims,id_RUG),fnam)
239  call nf95_def_var(nid, "longitude", NF90_FLOAT, ndim, varid_longitude)
240  call nf95_def_var(nid, "latitude",  NF90_FLOAT, ndim, varid_latitude)
241
242  !--- Attributes creation
243  CALL ncerr(NF90_PUT_ATT(nid,id_tim, "title","Jour dans l annee"),fnam)
244  CALL ncerr(NF90_PUT_ATT(nid,id_FOCE,"title","Fraction ocean"),fnam)
245  CALL ncerr(NF90_PUT_ATT(nid,id_FSIC,"title","Fraction glace de mer"),fnam)
246  CALL ncerr(NF90_PUT_ATT(nid,id_FTER,"title","Fraction terre"),fnam)
247  CALL ncerr(NF90_PUT_ATT(nid,id_FLIC,"title","Fraction land ice"),fnam)
248  CALL ncerr(NF90_PUT_ATT(nid,id_SST ,"title","Temperature superficielle de la mer"),fnam)
249  CALL ncerr(NF90_PUT_ATT(nid,id_BILS,"title","Reference flux de chaleur au sol"),fnam)
250  CALL ncerr(NF90_PUT_ATT(nid,id_ALB, "title","Albedo a la surface"),fnam)
251  CALL ncerr(NF90_PUT_ATT(nid,id_RUG, "title","Rugosite"),fnam)
252
253  call nf95_put_att(nid, varid_longitude, "standard_name", "longitude")
254  call nf95_put_att(nid, varid_longitude, "units", "degrees_east")
255
256  call nf95_put_att(nid, varid_latitude, "standard_name", "latitude")
257  call nf95_put_att(nid, varid_latitude, "units", "degrees_north")
258
259  CALL ncerr(NF90_ENDDEF(nid),fnam)
260
261  !--- Variables saving
262  CALL ncerr(NF90_PUT_VAR(nid,id_tim,[(REAL(k),k=1,ndays)]),fnam)
263  CALL ncerr(NF90_PUT_VAR(nid,id_FOCE,pctsrf_t(:,is_oce,:),[1,1],[klon,ndays]),fnam)
264  CALL ncerr(NF90_PUT_VAR(nid,id_FSIC,pctsrf_t(:,is_sic,:),[1,1],[klon,ndays]),fnam)
265  CALL ncerr(NF90_PUT_VAR(nid,id_FTER,pctsrf_t(:,is_ter,:),[1,1],[klon,ndays]),fnam)
266  CALL ncerr(NF90_PUT_VAR(nid,id_FLIC,pctsrf_t(:,is_lic,:),[1,1],[klon,ndays]),fnam)
267  CALL ncerr(NF90_PUT_VAR(nid,id_SST ,phy_sst(:,:),[1,1],[klon,ndays]),fnam)
268  CALL ncerr(NF90_PUT_VAR(nid,id_BILS,phy_bil(:,:),[1,1],[klon,ndays]),fnam)
269  CALL ncerr(NF90_PUT_VAR(nid,id_ALB ,phy_alb(:,:),[1,1],[klon,ndays]),fnam)
270  CALL ncerr(NF90_PUT_VAR(nid,id_RUG ,phy_rug(:,:),[1,1],[klon,ndays]),fnam)
271  call nf95_put_var(nid, varid_longitude, longitude_deg)
272  call nf95_put_var(nid, varid_latitude, latitude_deg)
273
274  CALL ncerr(NF90_CLOSE(nid),fnam)
275
276  CALL msg(6,'Ecriture du fichier limit : fin')
277
278  DEALLOCATE(pctsrf_t,phy_sst,phy_bil,phy_alb,phy_rug)
279
280
281!===============================================================================
282!
283  CONTAINS
284!
285!===============================================================================
286
287
288!-------------------------------------------------------------------------------
289!
290SUBROUTINE get_2Dfield(fnam, varname, mode, ndays, champo, flag, mask)
291!
292!-----------------------------------------------------------------------------
293! Comments:
294!   There are two assumptions concerning the NetCDF files, that are satisfied
295!   with files that are conforming NC convention:
296!     1) The last dimension of the variables used is the time record.
297!     2) Dimensional variables have the same names as corresponding dimensions.
298!-----------------------------------------------------------------------------
299  USE netcdf, ONLY: NF90_OPEN, NF90_INQ_VARID, NF90_INQUIRE_VARIABLE, &
300       NF90_CLOSE, NF90_INQ_DIMID, NF90_INQUIRE_DIMENSION, NF90_GET_VAR, &
301       NF90_GET_ATT
302  USE pchsp_95_m, only: pchsp_95
303  USE pchfe_95_m, only: pchfe_95
304  USE arth_m, only: arth
305  USE indice_sol_mod
306
307  IMPLICIT NONE
308  include "dimensions.h"
309  include "paramet.h"
310  include "comgeom2.h"
311!-----------------------------------------------------------------------------
312! Arguments:
313  CHARACTER(LEN=*),  INTENT(IN)     :: fnam     ! NetCDF file name
314  CHARACTER(LEN=10), INTENT(IN)     :: varname  ! NetCDF variable name
315  CHARACTER(LEN=3),  INTENT(IN)     :: mode     ! RUG, SIC, SST or ALB
316  INTEGER,           INTENT(IN)     :: ndays    ! current year number of days
317  REAL,    POINTER,  DIMENSION(:, :) :: champo  ! output field = f(t)
318  LOGICAL, OPTIONAL, INTENT(IN)     :: flag     ! extrapol. (SST) old ice (SIC)
319  REAL,    OPTIONAL, DIMENSION(iim, jjp1), INTENT(IN) :: mask
320!------------------------------------------------------------------------------
321! Local variables:
322!--- NetCDF
323  INTEGER           :: ncid, varid        ! NetCDF identifiers
324  CHARACTER(LEN=30) :: dnam               ! dimension name
325!--- dimensions
326  INTEGER           :: dids(4)            ! NetCDF dimensions identifiers
327  REAL, ALLOCATABLE :: dlon_ini(:)        ! initial longitudes vector
328  REAL, ALLOCATABLE :: dlat_ini(:)        ! initial latitudes  vector
329  REAL, POINTER     :: dlon(:), dlat(:)   ! reordered lon/lat  vectors
330!--- fields
331  INTEGER :: imdep, jmdep, lmdep          ! dimensions of 'champ'
332  REAL, ALLOCATABLE :: champ(:,:)         ! wanted field on initial grid
333  REAL, ALLOCATABLE :: yder(:), timeyear(:)
334  REAL              :: champint(iim,jjp1) ! interpolated field
335  REAL, ALLOCATABLE :: champtime(:,:,:)
336  REAL, ALLOCATABLE :: champan(:,:,:)
337!--- input files
338  CHARACTER(LEN=20) :: cal_in             ! calendar
339  CHARACTER(LEN=20) :: unit_sic           ! attribute unit in sea-ice file
340  INTEGER           :: ndays_in           ! number of days
341!--- misc
342  INTEGER           :: i, j, k, l, ll     ! loop counters
343  REAL, ALLOCATABLE :: work(:,:)          ! used for extrapolation
344  CHARACTER(LEN=128):: title, mess        ! for messages
345  LOGICAL           :: extrp              ! flag for extrapolation
346  REAL              :: chmin, chmax
347  INTEGER ierr
348  integer n_extrap ! number of extrapolated points
349  logical skip
350
351!------------------------------------------------------------------------------
352!---Variables depending on keyword 'mode' -------------------------------------
353  NULLIFY(champo)
354
355  SELECT CASE(mode)
356  CASE('RUG'); title='Rugosite'
357  CASE('SIC'); title='Sea-ice'
358  CASE('SST'); title='SST'
359  CASE('ALB'); title='Albedo'
360  END SELECT
361  extrp=.FALSE.; IF(PRESENT(flag).AND.mode=='SST') extrp=flag
362
363!--- GETTING SOME DIMENSIONAL VARIABLES FROM FILE -----------------------------
364  CALL msg(5,' Now reading file : '//TRIM(fnam))
365  CALL ncerr(NF90_OPEN(fnam, NF90_NOWRITE, ncid),fnam)
366  CALL ncerr(NF90_INQ_VARID(ncid, trim(varname), varid),fnam)
367  CALL ncerr(NF90_INQUIRE_VARIABLE(ncid, varid, dimids=dids),fnam)
368
369!--- Read unit for sea-ice variable only
370  IF (mode=='SIC') THEN
371    IF(NF90_GET_ATT(ncid, varid, 'units', unit_sic)/=NF90_NOERR) THEN
372      CALL msg(5,'No unit in sea-ice file. Take percentage as default value')
373      unit_sic='X'
374    ELSE
375      CALL msg(5,'Sea-ice cover has unit='//TRIM(unit_sic))
376    END IF
377  END IF
378
379!--- Longitude
380  CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(1), name=dnam, len=imdep),fnam)
381  ALLOCATE(dlon_ini(imdep), dlon(imdep))
382  CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam)
383  CALL ncerr(NF90_GET_VAR(ncid, varid, dlon_ini), fnam)
384  CALL msg(5,'variable '//TRIM(dnam)//' dimension ', imdep)
385
386!--- Latitude
387  CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(2), name=dnam, len=jmdep),fnam)
388  ALLOCATE(dlat_ini(jmdep), dlat(jmdep))
389  CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam)
390  CALL ncerr(NF90_GET_VAR(ncid, varid, dlat_ini), fnam)
391  CALL msg(5,'variable '//TRIM(dnam)//' dimension ', jmdep)
392
393!--- Time (variable is not needed - it is rebuilt - but calendar is)
394  CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(3), name=dnam, len=lmdep), fnam)
395  ALLOCATE(timeyear(MAX(lmdep,14)))
396  CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam)
397  cal_in=' '
398  IF(NF90_GET_ATT(ncid, varid, 'calendar', cal_in)/=NF90_NOERR) THEN
399    SELECT CASE(mode)
400      CASE('RUG', 'ALB'); cal_in='360d'
401      CASE('SIC', 'SST'); cal_in='gregorian'
402    END SELECT
403  CALL msg(5,'WARNING: missing "calendar" attribute for "time" in '&
404     &//TRIM(fnam)//'. Choosing default value.')
405  END IF
406  CALL msg(5,'var, calendar, dim: '//TRIM(dnam)//' '//TRIM(cal_in), lmdep)
407 
408!--- Determining input file number of days, depending on calendar
409  ndays_in=year_len(anneeref, cal_in)
410
411!--- Rebuilding input time vector (field from input file might be unreliable)
412  IF(lmdep==12) THEN
413    timeyear=mid_month(anneeref, cal_in, ndays_in)
414    CALL msg(0,'Monthly periodic input file (perpetual run).')
415  ELSE IF(lmdep==14) THEN
416    timeyear=mid_month(anneeref, cal_in, ndays_in)
417    CALL msg(0,'Monthly 14-records input file (interannual run).')
418  ELSE IF(lmdep==ndays_in) THEN
419    timeyear=[(REAL(k)-0.5,k=1,ndays_in)]
420    CALL msg(0,'Daily input file (no time interpolation).')
421  ELSE
422    WRITE(mess,'(a,i3,a,i3,a)')'Mismatching input file: found',lmdep,        &
423      ' records, 12/14/',ndays_in,' (periodic/interannual/daily) needed'
424    CALL abort_physic('mid_months',TRIM(mess),1)
425  END IF
426
427!--- GETTING THE FIELD AND INTERPOLATING IT ----------------------------------
428  ALLOCATE(champ(imdep, jmdep), champtime(iim, jjp1, MAX(lmdep,14)))
429  IF(extrp) ALLOCATE(work(imdep, jmdep))
430  CALL msg(5,'')
431  CALL msg(5,'LECTURE ET INTERPOLATION HORIZ. DE ', lmdep, ' CHAMPS.')
432  CALL ncerr(NF90_INQ_VARID(ncid, varname, varid), fnam)
433  DO l=1, lmdep
434    CALL ncerr(NF90_GET_VAR(ncid,varid,champ,[1,1,l],[imdep,jmdep,1]),fnam)
435    CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.)
436    IF(extrp) CALL extrapol(champ,imdep,jmdep,999999.,.TRUE.,.TRUE.,2,work)
437    IF(l==1) THEN
438      CALL msg(5,"----------------------------------------------------------")
439      CALL msg(5,"$$$ Interpolation barycentrique pour "//TRIM(title)//" $$$")
440      CALL msg(5,"----------------------------------------------------------")
441    END IF
442    IF(mode=='RUG') champ=LOG(champ)
443    CALL inter_barxy(dlon,dlat(:jmdep-1),champ,rlonu(:iim),rlatv,champint)
444    IF(mode=='RUG') THEN
445      champint=EXP(champint)
446      WHERE(NINT(mask)/=1) champint=0.001
447    END IF
448    !--- Special case for periodic input file: index shifted
449    ll = l; IF(lmdep==12) ll = l + 1
450    champtime(:, :, ll)=champint
451  END DO
452  IF(lmdep==12) THEN
453    champtime(:,:, 1)=champtime(:,:,13)
454    champtime(:,:,14)=champtime(:,:, 2)
455  END IF
456  CALL ncerr(NF90_CLOSE(ncid), fnam)
457
458  DEALLOCATE(dlon_ini, dlat_ini, dlon, dlat, champ)
459  IF(extrp) DEALLOCATE(work)
460
461!--- TIME INTERPOLATION ------------------------------------------------------
462  IF(prt_level>0) THEN
463     IF(ndays/=ndays_in) THEN
464        WRITE(lunout, *)
465        WRITE(lunout,*)'DIFFERENTES LONGEURS D ANNEES:'
466        WRITE(lunout,*)' Dans le fichier d entree: ',ndays_in
467        WRITE(lunout,*)' Dans le fichier de sortie: ',ndays
468     END IF
469     IF(lmdep==ndays_in) THEN
470        WRITE(lunout, *)
471        WRITE(lunout, *)'PAS D INTERPOLATION TEMPORELLE.'
472        WRITE(lunout, *)' Fichier journalier en entree.'
473     ELSE
474        WRITE(lunout, *)'INTERPOLATION TEMPORELLE.'
475        WRITE(lunout, *)' Vecteur temps en entree: ', timeyear
476        WRITE(lunout, *)' Vecteur temps en sortie de 0 a ', ndays-1
477     END IF
478  END IF
479  ALLOCATE(yder(MAX(lmdep,14)), champan(iip1, jjp1, ndays))
480  IF(lmdep==ndays_in) THEN
481     champan(1:iim,:,:)=champtime
482  ELSE
483     skip = .false.
484     n_extrap = 0
485     DO j=1, jjp1
486       DO i=1, iim
487         yder = pchsp_95(timeyear, champtime(i, j, :), ibeg=2, iend=2, &
488              vc_beg=0., vc_end=0.)
489         CALL pchfe_95(timeyear, champtime(i, j, :), yder, skip, &
490              arth(0., real(ndays_in) / ndays, ndays), champan(i, j, :), ierr)
491         if (ierr < 0) stop 1
492         n_extrap = n_extrap + ierr
493       END DO
494     END DO
495     IF(n_extrap /= 0) WRITE(lunout,*) "get_2Dfield pchfe_95: n_extrap = ", n_extrap
496  END IF
497  champan(iip1, :, :)=champan(1, :, :)
498  DEALLOCATE(yder, champtime, timeyear)
499
500!--- Checking the result
501  DO j=1, jjp1
502    CALL minmax(iip1, champan(1, j, 10), chmin, chmax)
503    IF (prt_level>5) WRITE(lunout, *)' ',TRIM(title),' au temps 10 ', chmin, chmax, j
504  END DO
505
506!--- SPECIAL FILTER FOR SST: SST>271.38 --------------------------------------
507  IF(mode=='SST') THEN
508    CALL msg(5,'Filtrage de la SST: SST >= 271.38')
509    WHERE(champan<271.38) champan=271.38
510  END IF
511
512!--- SPECIAL FILTER FOR SIC: 0.0<SIC<1.0 -------------------------------------
513  IF(mode=='SIC') THEN
514    CALL msg(5,'Filtrage de la SIC: 0.0 < Sea-ice < 1.0')
515
516    IF (unit_sic=='1') THEN
517       ! Nothing to be done for sea-ice field is already in fraction of 1
518       ! This is the case for sea-ice in file cpl_atm_sic.nc
519       CALL msg(5,'Sea-ice field already in fraction of 1')
520    ELSE
521       ! Convert sea ice from percentage to fraction of 1
522       CALL msg(5,'Transformt sea-ice field from percentage to fraction of 1.')
523       champan(:, :, :)=champan(:, :, :)/100.
524    END IF
525
526    champan(iip1, :, :)=champan(1, :, :)
527    WHERE(champan>1.0) champan=1.0
528    WHERE(champan<0.0) champan=0.0
529 END IF
530
531!--- DYNAMICAL TO PHYSICAL GRID ----------------------------------------------
532  ALLOCATE(champo(klon, ndays))
533  DO k=1, ndays
534    CALL gr_dyn_fi(1, iip1, jjp1, klon, champan(1, 1, k), champo(1, k))
535  END DO
536  DEALLOCATE(champan)
537
538END SUBROUTINE get_2Dfield
539!
540!-------------------------------------------------------------------------------
541
542
543!-------------------------------------------------------------------------------
544!
545SUBROUTINE start_init_orog0(lon_in,lat_in,phis,masque)
546!
547!-------------------------------------------------------------------------------
548  IMPLICIT NONE
549!===============================================================================
550! Purpose:  Compute "phis" just like it would be in start_init_orog.
551!===============================================================================
552! Arguments:
553  REAL,             INTENT(IN)    :: lon_in(:), lat_in(:)   ! dim (iml) (jml)
554  REAL,             INTENT(INOUT) :: phis(:,:), masque(:,:) ! dim (iml,jml)
555!-------------------------------------------------------------------------------
556! Local variables:
557  CHARACTER(LEN=256) :: modname="start_init_orog0"
558  INTEGER            :: fid, llm_tmp,ttm_tmp, iml,jml, iml_rel,jml_rel, itau(1)
559  REAL               :: lev(1), date, dt, deg2rad
560  REAL, ALLOCATABLE  :: lon_rad(:), lon_ini(:), lon_rel(:,:), relief_hi(:,:)
561  REAL, ALLOCATABLE  :: lat_rad(:), lat_ini(:), lat_rel(:,:)
562!-------------------------------------------------------------------------------
563  iml=assert_eq(SIZE(lon_in),SIZE(phis,1),SIZE(masque,1),TRIM(modname)//" iml")
564  jml=assert_eq(SIZE(lat_in),SIZE(phis,2),SIZE(masque,2),TRIM(modname)//" jml")
565  IF(iml/=iip1) CALL abort_gcm(TRIM(modname),'iml/=iip1',1)
566  IF(jml/=jjp1) CALL abort_gcm(TRIM(modname),'jml/=jjp1',1)
567  pi=2.0*ASIN(1.0); deg2rad=pi/180.0
568  IF(ANY(phis/=-99999.)) RETURN                  !--- phis ALREADY KNOWN
569
570!--- HIGH RESOLUTION OROGRAPHY
571  CALL flininfo(frelf, iml_rel, jml_rel, llm_tmp, ttm_tmp, fid)
572
573  ALLOCATE(lat_rel(iml_rel,jml_rel),lon_rel(iml_rel,jml_rel))
574  CALL flinopen(frelf, .FALSE., iml_rel, jml_rel, llm_tmp, lon_rel, lat_rel,   &
575                lev, ttm_tmp, itau, date, dt, fid)
576  ALLOCATE(relief_hi(iml_rel,jml_rel))
577  CALL flinget(fid, vrel, iml_rel, jml_rel, llm_tmp, ttm_tmp, 1, 1, relief_hi)
578  CALL flinclo(fid)
579
580!--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS
581  ALLOCATE(lon_ini(iml_rel),lat_ini(jml_rel))
582  lon_ini(:)=lon_rel(:,1); IF(MAXVAL(lon_rel)>pi) lon_ini=lon_ini*deg2rad
583  lat_ini(:)=lat_rel(1,:); IF(MAXVAL(lat_rel)>pi) lat_ini=lat_ini*deg2rad
584
585!--- FIELDS ARE PROCESSED TO BE ON STANDARD ANGULAR DOMAINS
586  ALLOCATE(lon_rad(iml_rel),lat_rad(jml_rel))
587  CALL conf_dat2d(vrel, lon_ini, lat_ini, lon_rad, lat_rad, relief_hi, .FALSE.)
588  DEALLOCATE(lon_ini,lat_ini)
589
590!--- COMPUTING SURFACE GEOPOTENTIAL USING ROUTINE grid_noro0
591  WRITE(lunout,*)
592  WRITE(lunout,*)'*** Compute surface geopotential ***'
593
594!--- CALL OROGRAPHY MODULE (REDUCED VERSION) TO COMPUTE FIELDS
595  CALL grid_noro0(lon_rad, lat_rad, relief_hi, lon_in, lat_in, phis, masque)
596  phis = phis * 9.81
597  phis(iml,:) = phis(1,:)
598  DEALLOCATE(relief_hi,lon_rad,lat_rad)
599
600END SUBROUTINE start_init_orog0
601!
602!-------------------------------------------------------------------------------
603
604
605!-------------------------------------------------------------------------------
606!
607SUBROUTINE grid_noro0(xd,yd,zd,x,y,zphi,mask)
608!
609!===============================================================================
610! Purpose: Extracted from grid_noro to provide geopotential height for dynamics
611!          without any call to physics subroutines.
612!===============================================================================
613  IMPLICIT NONE
614!-------------------------------------------------------------------------------
615! Arguments:
616  REAL, INTENT(IN)   :: xd(:), yd(:) !--- INPUT  COORDINATES     (imdp) (jmdp)
617  REAL, INTENT(IN)   :: zd(:,:)      !--- INPUT  FIELD           (imdp,jmdp)
618  REAL, INTENT(IN)   :: x(:), y(:)   !--- OUTPUT COORDINATES     (imar+1) (jmar)
619  REAL, INTENT(OUT)  :: zphi(:,:)    !--- GEOPOTENTIAL           (imar+1,jmar)
620  REAL, INTENT(INOUT):: mask(:,:)    !--- MASK                   (imar+1,jmar)
621!-------------------------------------------------------------------------------
622! Local variables:
623  CHARACTER(LEN=256) :: modname="grid_noro0"
624  REAL, ALLOCATABLE :: xusn(:), yusn(:)           ! dim (imdp+2*iext) (jmdp+2)
625  REAL, ALLOCATABLE :: zusn(:,:)                  ! dim (imdp+2*iext,jmdp+2)
626  REAL, ALLOCATABLE :: weight(:,:)                ! dim (imar+1,jmar)
627  REAL, ALLOCATABLE :: mask_tmp(:,:), zmea(:,:)   ! dim (imar+1,jmar)
628  REAL, ALLOCATABLE :: num_tot(:,:), num_lan(:,:) ! dim (imax,jmax)
629  REAL, ALLOCATABLE :: a(:), b(:)                 ! dim (imax)
630  REAL, ALLOCATABLE :: c(:), d(:)                 ! dim (jmax)
631  LOGICAL :: masque_lu
632  INTEGER :: i, ii, imdp, imar, iext
633  INTEGER :: j, jj, jmdp, jmar, nn
634  REAL    :: xpi, zlenx, weighx, xincr,  zbordnor, zmeanor, zweinor, zbordest
635  REAL    :: rad, zleny, weighy, masque, zbordsud, zmeasud, zweisud, zbordoue
636!-------------------------------------------------------------------------------
637  imdp=assert_eq(SIZE(xd),SIZE(zd,1),TRIM(modname)//" imdp")
638  jmdp=assert_eq(SIZE(yd),SIZE(zd,2),TRIM(modname)//" jmdp")
639  imar=assert_eq(SIZE(x),SIZE(zphi,1),SIZE(mask,1),TRIM(modname)//" imar")-1
640  jmar=assert_eq(SIZE(y),SIZE(zphi,2),SIZE(mask,2),TRIM(modname)//" jmar")
641  IF(imar/=iim)   CALL abort_gcm(TRIM(modname),'imar/=iim'  ,1)
642  IF(jmar/=jjm+1) CALL abort_gcm(TRIM(modname),'jmar/=jjm+1',1)
643  iext=imdp/10
644  xpi = ACOS(-1.)
645  rad = 6371229.
646
647!--- ARE WE USING A READ MASK ?
648  masque_lu=ANY(mask/=-99999.); IF(.NOT.masque_lu) mask=0.0
649  WRITE(lunout,*)'Masque lu: ',masque_lu
650
651!--- EXTENSION OF THE INPUT DATABASE TO PROCEED COMPUTATIONS AT BOUNDARIES:
652  ALLOCATE(xusn(imdp+2*iext))
653  xusn(1     +iext:imdp  +iext)=xd(:)
654  xusn(1          :       iext)=xd(1+imdp-iext:imdp)-2.*xpi
655  xusn(1+imdp+iext:imdp+2*iext)=xd(1          :iext)+2.*xpi
656
657  ALLOCATE(yusn(jmdp+2))
658  yusn(1       )=yd(1)   +(yd(1)   -yd(2))
659  yusn(2:jmdp+1)=yd(:)
660  yusn(  jmdp+2)=yd(jmdp)+(yd(jmdp)-yd(jmdp-1))
661
662  ALLOCATE(zusn(imdp+2*iext,jmdp+2))
663  zusn(1       +iext:imdp  +iext,2:jmdp+1)=zd  (:                   ,     :)
664  zusn(1            :       iext,2:jmdp+1)=zd  (imdp-iext+1:imdp    ,     :)
665  zusn(1+imdp  +iext:imdp+2*iext,2:jmdp+1)=zd  (1:iext              ,     :)
666  zusn(1            :imdp/2+iext,       1)=zusn(1+imdp/2:imdp  +iext,     2)
667  zusn(1+imdp/2+iext:imdp+2*iext,       1)=zusn(1       :imdp/2+iext,     2)
668  zusn(1            :imdp/2+iext,  jmdp+2)=zusn(1+imdp/2:imdp  +iext,jmdp+1)
669  zusn(1+imdp/2+iext:imdp+2*iext,  jmdp+2)=zusn(1       :imdp/2+iext,jmdp+1)
670
671!--- COMPUTE LIMITS OF MODEL GRIDPOINT AREA (REGULAR GRID)
672  ALLOCATE(a(imar+1),b(imar+1))
673  b(1:imar)=(x(1:imar  )+ x(2:imar+1))/2.0
674  b(imar+1)= x(  imar+1)+(x(  imar+1)-x(imar))/2.0
675  a(1)=x(1)-(x(2)-x(1))/2.0
676  a(2:imar+1)= b(1:imar)
677
678  ALLOCATE(c(jmar),d(jmar))
679  d(1:jmar-1)=(y(1:jmar-1)+ y(2:jmar))/2.0
680  d(  jmar  )= y(  jmar  )+(y(  jmar)-y(jmar-1))/2.0
681  c(1)=y(1)-(y(2)-y(1))/2.0
682  c(2:jmar)=d(1:jmar-1)
683
684!--- INITIALIZATIONS:
685  ALLOCATE(weight(imar+1,jmar)); weight(:,:)= 0.0
686  ALLOCATE(zmea  (imar+1,jmar)); zmea  (:,:)= 0.0
687
688!--- SUMMATION OVER GRIDPOINT AREA
689  zleny=xpi/REAL(jmdp)*rad
690  xincr=xpi/REAL(jmdp)/2.
691  ALLOCATE(num_tot(imar+1,jmar)); num_tot(:,:)=0.
692  ALLOCATE(num_lan(imar+1,jmar)); num_lan(:,:)=0.
693  DO ii = 1, imar+1
694    DO jj = 1, jmar
695      DO j = 2,jmdp+1
696        zlenx  =zleny  *COS(yusn(j))
697        zbordnor=(xincr+c(jj)-yusn(j))*rad
698        zbordsud=(xincr-d(jj)+yusn(j))*rad
699        weighy=AMAX1(0.,AMIN1(zbordnor,zbordsud,zleny))
700        IF(weighy/=0) THEN
701          DO i = 2, imdp+2*iext-1
702            zbordest=(xusn(i)-a(ii)+xincr)*rad*COS(yusn(j))
703            zbordoue=(b(ii)+xincr-xusn(i))*rad*COS(yusn(j))
704            weighx=AMAX1(0.,AMIN1(zbordest,zbordoue,zlenx))
705            IF(weighx/=0)THEN
706              num_tot(ii,jj)=num_tot(ii,jj)+1.0
707              IF(zusn(i,j)>=1.)num_lan(ii,jj)=num_lan(ii,jj)+1.0
708              weight(ii,jj)=weight(ii,jj)+weighx*weighy
709              zmea  (ii,jj)=zmea  (ii,jj)+zusn(i,j)*weighx*weighy !--- MEAN
710            END IF
711          END DO
712        END IF
713      END DO
714    END DO
715  END DO
716
717!--- COMPUTE PARAMETERS NEEDED BY LOTT & MILLER (1997) AND LOTT (1999) SSO SCHEME
718  IF(.NOT.masque_lu) THEN
719    WHERE(weight(:,1:jmar-1)/=0.0) mask=num_lan(:,:)/num_tot(:,:)
720  END IF
721  nn=COUNT(weight(:,1:jmar-1)==0.0)
722  IF(nn/=0) WRITE(lunout,*)'Problem with weight ; vanishing occurrences: ',nn
723  WHERE(weight/=0.0) zmea(:,:)=zmea(:,:)/weight(:,:)
724
725!--- MASK BASED ON GROUND MAXIMUM, 10% THRESHOLD (<10%: SURF PARAMS MEANINGLESS)
726  ALLOCATE(mask_tmp(imar+1,jmar)); mask_tmp(:,:)=0.0
727  WHERE(mask>=0.1) mask_tmp = 1.
728  WHERE(weight(:,:)/=0.0)
729    zphi(:,:)=mask_tmp(:,:)*zmea(:,:)
730    zmea(:,:)=mask_tmp(:,:)*zmea(:,:)
731  END WHERE
732  WRITE(lunout,*)'  MEAN ORO:' ,MAXVAL(zmea)
733
734!--- Values at poles
735  zphi(imar+1,:)=zphi(1,:)
736
737  zweinor=SUM(weight(1:imar,   1),DIM=1)
738  zweisud=SUM(weight(1:imar,jmar),DIM=1)
739  zmeanor=SUM(weight(1:imar,   1)*zmea(1:imar,   1),DIM=1)
740  zmeasud=SUM(weight(1:imar,jmar)*zmea(1:imar,jmar),DIM=1)
741  zphi(:,1)=zmeanor/zweinor; zphi(:,jmar)=zmeasud/zweisud
742
743END SUBROUTINE grid_noro0
744!
745!-------------------------------------------------------------------------------
746
747
748!-------------------------------------------------------------------------------
749!
750FUNCTION year_len(y,cal_in)
751!
752!-------------------------------------------------------------------------------
753  IMPLICIT NONE
754!-------------------------------------------------------------------------------
755! Arguments:
756  INTEGER                       :: year_len
757  INTEGER,           INTENT(IN) :: y
758  CHARACTER(LEN=*),  INTENT(IN) :: cal_in
759!-------------------------------------------------------------------------------
760! Local variables:
761  CHARACTER(LEN=20)             :: cal_out              ! calendar (for outputs)
762!-------------------------------------------------------------------------------
763!--- Getting the input calendar to reset at the end of the function
764  CALL ioget_calendar(cal_out)
765
766!--- Unlocking calendar and setting it to wanted one
767  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_in))
768
769!--- Getting the number of days in this year
770  year_len=ioget_year_len(y)
771
772!--- Back to original calendar
773  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_out))
774
775END FUNCTION year_len
776!
777!-------------------------------------------------------------------------------
778
779
780!-------------------------------------------------------------------------------
781!
782FUNCTION mid_month(y,cal_in,ndays_out)
783!
784!-------------------------------------------------------------------------------
785  IMPLICIT NONE
786!-------------------------------------------------------------------------------
787! Arguments:
788  INTEGER,          INTENT(IN) :: y           ! year
789  CHARACTER(LEN=*), INTENT(IN) :: cal_in      ! calendar
790  INTEGER,          INTENT(IN) :: ndays_out   ! days   number
791  REAL,          DIMENSION(14) :: mid_month   ! mid-bins times
792!-------------------------------------------------------------------------------
793! Local variables:
794  CHARACTER(LEN=99)      :: mess              ! error message
795  CHARACTER(LEN=20)      :: cal_out           ! output calendar
796  INTEGER, DIMENSION(14) :: tlen              ! months lengths (days)
797  INTEGER                :: m                 ! months counter
798  INTEGER                :: nd                ! number of days
799!-------------------------------------------------------------------------------
800  !--- Save the input calendar to reset it at the end of the function
801  CALL ioget_calendar(cal_out)
802
803  !--- Unlock calendar and set it to "cal_in"
804  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_in))
805
806  !--- Get the length of each month
807  tlen(1 )=ioget_mon_len(y-1,12)
808  DO m=1,12; tlen(m+1)=ioget_mon_len(y,m); END DO
809  tlen(14)=ioget_mon_len(y+1, 1)
810
811  !--- Mid-bins times
812  mid_month(1)=-0.5*REAL(tlen(1))
813  DO m=2,14; mid_month(m)=mid_month(m-1)+0.5*REAL(tlen(m-1)+tlen(m)); END DO
814
815  !--- Back to original calendar
816  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_out))
817
818END FUNCTION mid_month
819!
820!-------------------------------------------------------------------------------
821
822
823
824!-------------------------------------------------------------------------------
825!
826SUBROUTINE msg(lev,str1,i,str2)
827!
828!-------------------------------------------------------------------------------
829! Arguments:
830  INTEGER,                    INTENT(IN) :: lev
831  CHARACTER(LEN=*),           INTENT(IN) :: str1
832  INTEGER,          OPTIONAL, INTENT(IN) :: i
833  CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: str2
834!-------------------------------------------------------------------------------
835  IF(prt_level>lev) THEN
836    IF(PRESENT(str2)) THEN
837      WRITE(lunout,*) TRIM(str1), i, TRIM(str2)
838    ELSE IF(PRESENT(i)) THEN
839      WRITE(lunout,*) TRIM(str1), i
840    ELSE
841      WRITE(lunout,*) TRIM(str1)
842    END IF
843  END IF
844
845END SUBROUTINE msg
846!
847!-------------------------------------------------------------------------------
848
849
850!-------------------------------------------------------------------------------
851!
852SUBROUTINE ncerr(ncres,fnam)
853!
854!-------------------------------------------------------------------------------
855! Purpose: NetCDF errors handling.
856!-------------------------------------------------------------------------------
857  USE netcdf, ONLY : NF90_NOERR, NF90_STRERROR
858  IMPLICIT NONE
859!-------------------------------------------------------------------------------
860! Arguments:
861  INTEGER,          INTENT(IN) :: ncres
862  CHARACTER(LEN=*), INTENT(IN) :: fnam
863!-------------------------------------------------------------------------------
864  IF(ncres/=NF90_NOERR) THEN
865    WRITE(lunout,*)'Problem with file '//TRIM(fnam)//' in routine limit_netcdf.'
866    CALL abort_physic('limit_netcdf',NF90_STRERROR(ncres),1)
867  END IF
868
869END SUBROUTINE ncerr
870!
871!-------------------------------------------------------------------------------
872
873#endif
874! of #ifndef CPP_1D
875END SUBROUTINE limit_netcdf
876
877END MODULE limit
878!
879!*******************************************************************************
880
Note: See TracBrowser for help on using the repository browser.