source: LMDZ5/trunk/libf/dynlonlat_phylonlat/phylmd/limit_netcdf.F90 @ 2342

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

Follow-up for a cleaner separation between dynamics and physics:

  • "write_field" is called from physics and dynamics but has no dependence on either so it should be in "misc"
  • "write_field_phy" is common to all physics, so it goes in "phy_common"
  • "init_ssrf_m" and "limit_netcdf" are only used by ce0l, so these should be in "dynlonlat_phylonlat/phylmd"
  • "q_sat" is called from physics and dynamics but has no dependence on either so it should be in "misc"

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