source: LMDZ5/branches/AI-cosp/libf/dynphy_lonlat/phylmd/limit_netcdf.F90 @ 3906

Last change on this file since 3906 was 2399, checked in by Ehouarn Millour, 9 years ago

Follow-up from commit 2395: get rid of rlon and rlat, longitude_deg and latitude_deg (from module geometry_mod) should be used instead. Longitudes and latitudes are no longer loaded from startphy.nc but inherited from dynamics (and compatibility with values in startphy.nc is checked). This will change bench results because of roundoffs differences between the two.
EM

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