source: LMDZ5/branches/IPSLCM5A2.1/libf/dynphy_lonlat/phylmd/limit_netcdf.F90 @ 5425

Last change on this file since 5425 was 2603, checked in by Ehouarn Millour, 8 years ago

Cleanup in the dynamics: turn logic.h into module logic_mod.F90
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: 32.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  USE comconst_mod, ONLY: pi
73  IMPLICIT NONE
74!-------------------------------------------------------------------------------
75! Arguments:
76  include "iniprint.h"
77  include "dimensions.h"
78  include "paramet.h"
79  REAL, DIMENSION(iip1,jjp1), INTENT(INOUT) :: masque ! land mask
80  REAL, DIMENSION(iip1,jjp1), INTENT(INOUT) :: phis   ! ground geopotential
81  LOGICAL,                    INTENT(IN)    :: extrap ! SST extrapolation flag
82!-------------------------------------------------------------------------------
83! Local variables:
84  include "comgeom2.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  !--- TO MATCH EXACTLY WHAT WOULD BE DONE IN etat0phys_netcdf
119    WHERE(   masque(:,:)<EPSFRA) masque(:,:)=0.
120    WHERE(1.-masque(:,:)<EPSFRA) masque(:,:)=1.
121  END IF
122
123!--- Beware: anneeref (from gcm.def) is used to determine output time sampling
124  ndays=ioget_year_len(anneeref)
125
126!--- RUGOSITY TREATMENT --------------------------------------------------------
127  CALL msg(1,'Traitement de la rugosite')
128  CALL get_2Dfield(frugo,vrug,'RUG',ndays,phy_rug,mask=masque(1:iim,:))
129
130!--- OCEAN TREATMENT -----------------------------------------------------------
131  CALL msg(1,'Traitement de la glace oceanique')
132
133! Input SIC file selection
134! Open file only to test if available
135  DO ix_sic=1,SIZE(fsic)
136     IF ( NF90_OPEN(TRIM(fsic(ix_sic)),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
137        icefile=fsic(ix_sic); varname=vsic(ix_sic); EXIT
138     END IF
139  END DO
140  IF(ix_sic==SIZE(fsic)+1) THEN
141     WRITE(lunout,*) 'ERROR! No sea-ice input file was found.'
142     WRITE(lunout,*) 'One of following files must be available : '
143     DO k=1,SIZE(fsic); WRITE(lunout,*) TRIM(fsic(k)); END DO
144     CALL abort_physic('limit_netcdf','No sea-ice file was found',1)
145  END IF
146  CALL ncerr(NF90_CLOSE(nid),icefile)
147  CALL msg(-1,'Fichier choisi pour la glace de mer:'//TRIM(icefile))
148
149  CALL get_2Dfield(icefile,varname, 'SIC',ndays,phy_ice)
150
151  ALLOCATE(pctsrf_t(klon,nbsrf,ndays))
152  DO k=1,ndays
153     fi_ice=phy_ice(:,k)
154     WHERE(fi_ice>=1.0  ) fi_ice=1.0
155     WHERE(fi_ice<EPSFRA) fi_ice=0.0
156     pctsrf_t(:,is_ter,k)=pctsrf(:,is_ter)       ! land soil
157     pctsrf_t(:,is_lic,k)=pctsrf(:,is_lic)       ! land ice
158     SELECT CASE(ix_sic)
159        CASE(2)                                   ! SIC=pICE*(1-LIC-TER) (CPL)
160        pctsrf_t(:,is_sic,k)=fi_ice(:)*(1.-pctsrf(:,is_lic)-pctsrf(:,is_ter))
161        CASE(3)                                   ! SIC=pICE            (HIST)
162        pctsrf_t(:,is_sic,k)=fi_ice(:)
163        CASE DEFAULT                              ! SIC=pICE-LIC   (AMIP,ERAI)
164        pctsrf_t(:,is_sic,k)=fi_ice-pctsrf_t(:,is_lic,k)
165     END SELECT
166     WHERE(pctsrf_t(:,is_sic,k)<=0) pctsrf_t(:,is_sic,k)=0.
167     WHERE(1.0-zmasq<EPSFRA)
168        pctsrf_t(:,is_sic,k)=0.0
169        pctsrf_t(:,is_oce,k)=0.0
170     ELSEWHERE
171        WHERE(pctsrf_t(:,is_sic,k)>=1.0-zmasq)
172           pctsrf_t(:,is_sic,k)=1.0-zmasq
173           pctsrf_t(:,is_oce,k)=0.0
174        ELSEWHERE
175           pctsrf_t(:,is_oce,k)=1.0-zmasq-pctsrf_t(:,is_sic,k)
176           WHERE(pctsrf_t(:,is_oce,k)<EPSFRA)
177              pctsrf_t(:,is_oce,k)=0.0
178              pctsrf_t(:,is_sic,k)=1.0-zmasq
179           END WHERE
180        END WHERE
181     END WHERE
182     nbad=COUNT(pctsrf_t(:,is_oce,k)<0.0)
183     IF(nbad>0) WRITE(lunout,*) 'pb sous maille pour nb points = ',nbad
184     nbad=COUNT(ABS(SUM(pctsrf_t(:,:,k),DIM=2)-1.0)>EPSFRA)
185     IF(nbad>0) WRITE(lunout,*) 'pb sous surface pour nb points = ',nbad
186  END DO
187  DEALLOCATE(phy_ice)
188
189!--- SST TREATMENT -------------------------------------------------------------
190  CALL msg(1,'Traitement de la sst')
191
192! Input SST file selection
193! Open file only to test if available
194  DO ix_sst=1,SIZE(fsst)
195     IF ( NF90_OPEN(TRIM(fsst(ix_sst)),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
196       sstfile=fsst(ix_sst); varname=vsst(ix_sst); EXIT
197     END IF
198  END DO
199  IF(ix_sst==SIZE(fsst)+1) THEN
200     WRITE(lunout,*) 'ERROR! No sst input file was found.'
201     WRITE(lunout,*) 'One of following files must be available : '
202     DO k=1,SIZE(fsst); WRITE(lunout,*) TRIM(fsst(k)); END DO
203     CALL abort_physic('limit_netcdf','No sst file was found',1)
204  END IF
205  CALL ncerr(NF90_CLOSE(nid),sstfile)
206  CALL msg(-1,'Fichier choisi pour la temperature de mer: '//TRIM(sstfile))
207
208  CALL get_2Dfield(sstfile,varname,'SST',ndays,phy_sst,flag=extrap)
209
210!--- ALBEDO TREATMENT ----------------------------------------------------------
211  CALL msg(1,'Traitement de l albedo')
212  CALL get_2Dfield(falbe,valb,'ALB',ndays,phy_alb)
213
214!--- REFERENCE GROUND HEAT FLUX TREATMENT --------------------------------------
215  ALLOCATE(phy_bil(klon,ndays)); phy_bil=0.0
216
217!--- OUTPUT FILE WRITING -------------------------------------------------------
218  CALL msg(5,'Ecriture du fichier limit : debut')
219  fnam="limit.nc"
220
221  !--- File creation
222  CALL ncerr(NF90_CREATE(fnam,NF90_CLOBBER,nid),fnam)
223  CALL ncerr(NF90_PUT_ATT(nid,NF90_GLOBAL,"title","Fichier conditions aux limites"),fnam)
224
225  !--- Dimensions creation
226  CALL ncerr(NF90_DEF_DIM(nid,"points_physiques",klon,ndim),fnam)
227  CALL ncerr(NF90_DEF_DIM(nid,"time",NF90_UNLIMITED,ntim),fnam)
228
229  dims=[ndim,ntim]
230
231  !--- Variables creation
232  CALL ncerr(NF90_DEF_VAR(nid,"TEMPS",NF90_FORMAT,[ntim],id_tim),fnam)
233  CALL ncerr(NF90_DEF_VAR(nid,"FOCE", NF90_FORMAT,dims,id_FOCE),fnam)
234  CALL ncerr(NF90_DEF_VAR(nid,"FSIC", NF90_FORMAT,dims,id_FSIC),fnam)
235  CALL ncerr(NF90_DEF_VAR(nid,"FTER", NF90_FORMAT,dims,id_FTER),fnam)
236  CALL ncerr(NF90_DEF_VAR(nid,"FLIC", NF90_FORMAT,dims,id_FLIC),fnam)
237  CALL ncerr(NF90_DEF_VAR(nid,"SST",  NF90_FORMAT,dims,id_SST),fnam)
238  CALL ncerr(NF90_DEF_VAR(nid,"BILS", NF90_FORMAT,dims,id_BILS),fnam)
239  CALL ncerr(NF90_DEF_VAR(nid,"ALB",  NF90_FORMAT,dims,id_ALB),fnam)
240  CALL ncerr(NF90_DEF_VAR(nid,"RUG",  NF90_FORMAT,dims,id_RUG),fnam)
241  call nf95_def_var(nid, "longitude", NF90_FLOAT, ndim, varid_longitude)
242  call nf95_def_var(nid, "latitude",  NF90_FLOAT, ndim, varid_latitude)
243
244  !--- Attributes creation
245  CALL ncerr(NF90_PUT_ATT(nid,id_tim, "title","Jour dans l annee"),fnam)
246  CALL ncerr(NF90_PUT_ATT(nid,id_FOCE,"title","Fraction ocean"),fnam)
247  CALL ncerr(NF90_PUT_ATT(nid,id_FSIC,"title","Fraction glace de mer"),fnam)
248  CALL ncerr(NF90_PUT_ATT(nid,id_FTER,"title","Fraction terre"),fnam)
249  CALL ncerr(NF90_PUT_ATT(nid,id_FLIC,"title","Fraction land ice"),fnam)
250  CALL ncerr(NF90_PUT_ATT(nid,id_SST ,"title","Temperature superficielle de la mer"),fnam)
251  CALL ncerr(NF90_PUT_ATT(nid,id_BILS,"title","Reference flux de chaleur au sol"),fnam)
252  CALL ncerr(NF90_PUT_ATT(nid,id_ALB, "title","Albedo a la surface"),fnam)
253  CALL ncerr(NF90_PUT_ATT(nid,id_RUG, "title","Rugosite"),fnam)
254
255  call nf95_put_att(nid, varid_longitude, "standard_name", "longitude")
256  call nf95_put_att(nid, varid_longitude, "units", "degrees_east")
257
258  call nf95_put_att(nid, varid_latitude, "standard_name", "latitude")
259  call nf95_put_att(nid, varid_latitude, "units", "degrees_north")
260
261  CALL ncerr(NF90_ENDDEF(nid),fnam)
262
263  !--- Variables saving
264  CALL ncerr(NF90_PUT_VAR(nid,id_tim,[(REAL(k),k=1,ndays)]),fnam)
265  CALL ncerr(NF90_PUT_VAR(nid,id_FOCE,pctsrf_t(:,is_oce,:),[1,1],[klon,ndays]),fnam)
266  CALL ncerr(NF90_PUT_VAR(nid,id_FSIC,pctsrf_t(:,is_sic,:),[1,1],[klon,ndays]),fnam)
267  CALL ncerr(NF90_PUT_VAR(nid,id_FTER,pctsrf_t(:,is_ter,:),[1,1],[klon,ndays]),fnam)
268  CALL ncerr(NF90_PUT_VAR(nid,id_FLIC,pctsrf_t(:,is_lic,:),[1,1],[klon,ndays]),fnam)
269  CALL ncerr(NF90_PUT_VAR(nid,id_SST ,phy_sst(:,:),[1,1],[klon,ndays]),fnam)
270  CALL ncerr(NF90_PUT_VAR(nid,id_BILS,phy_bil(:,:),[1,1],[klon,ndays]),fnam)
271  CALL ncerr(NF90_PUT_VAR(nid,id_ALB ,phy_alb(:,:),[1,1],[klon,ndays]),fnam)
272  CALL ncerr(NF90_PUT_VAR(nid,id_RUG ,phy_rug(:,:),[1,1],[klon,ndays]),fnam)
273  call nf95_put_var(nid, varid_longitude, longitude_deg)
274  call nf95_put_var(nid, varid_latitude, latitude_deg)
275
276  CALL ncerr(NF90_CLOSE(nid),fnam)
277
278  CALL msg(6,'Ecriture du fichier limit : fin')
279
280  DEALLOCATE(pctsrf_t,phy_sst,phy_bil,phy_alb,phy_rug)
281
282
283!===============================================================================
284!
285  CONTAINS
286!
287!===============================================================================
288
289
290!-------------------------------------------------------------------------------
291!
292SUBROUTINE get_2Dfield(fnam, varname, mode, ndays, champo, flag, mask)
293!
294!-----------------------------------------------------------------------------
295! Comments:
296!   There are two assumptions concerning the NetCDF files, that are satisfied
297!   with files that are conforming NC convention:
298!     1) The last dimension of the variables used is the time record.
299!     2) Dimensional variables have the same names as corresponding dimensions.
300!-----------------------------------------------------------------------------
301  USE netcdf, ONLY: NF90_OPEN, NF90_INQ_VARID, NF90_INQUIRE_VARIABLE, &
302       NF90_CLOSE, NF90_INQ_DIMID, NF90_INQUIRE_DIMENSION, NF90_GET_VAR, &
303       NF90_GET_ATT
304  USE pchsp_95_m, only: pchsp_95
305  USE pchfe_95_m, only: pchfe_95
306  USE arth_m, only: arth
307  USE indice_sol_mod
308
309  IMPLICIT NONE
310  include "dimensions.h"
311  include "paramet.h"
312  include "comgeom2.h"
313!-----------------------------------------------------------------------------
314! Arguments:
315  CHARACTER(LEN=*),  INTENT(IN)     :: fnam     ! NetCDF file name
316  CHARACTER(LEN=10), INTENT(IN)     :: varname  ! NetCDF variable name
317  CHARACTER(LEN=3),  INTENT(IN)     :: mode     ! RUG, SIC, SST or ALB
318  INTEGER,           INTENT(IN)     :: ndays    ! current year number of days
319  REAL,    POINTER,  DIMENSION(:, :) :: champo  ! output field = f(t)
320  LOGICAL, OPTIONAL, INTENT(IN)     :: flag     ! extrapol. (SST) old ice (SIC)
321  REAL,    OPTIONAL, DIMENSION(iim, jjp1), INTENT(IN) :: mask
322!------------------------------------------------------------------------------
323! Local variables:
324!--- NetCDF
325  INTEGER           :: ncid, varid        ! NetCDF identifiers
326  CHARACTER(LEN=30) :: dnam               ! dimension name
327!--- dimensions
328  INTEGER           :: dids(4)            ! NetCDF dimensions identifiers
329  REAL, ALLOCATABLE :: dlon_ini(:)        ! initial longitudes vector
330  REAL, ALLOCATABLE :: dlat_ini(:)        ! initial latitudes  vector
331  REAL, POINTER     :: dlon(:), dlat(:)   ! reordered lon/lat  vectors
332!--- fields
333  INTEGER :: imdep, jmdep, lmdep          ! dimensions of 'champ'
334  REAL, ALLOCATABLE :: champ(:,:)         ! wanted field on initial grid
335  REAL, ALLOCATABLE :: yder(:), timeyear(:)
336  REAL              :: champint(iim,jjp1) ! interpolated field
337  REAL, ALLOCATABLE :: champtime(:,:,:)
338  REAL, ALLOCATABLE :: champan(:,:,:)
339!--- input files
340  CHARACTER(LEN=20) :: fnam_m, fnam_p     ! previous/next files names
341  CHARACTER(LEN=20) :: cal_in             ! calendar
342  CHARACTER(LEN=20) :: unit_sic           ! attribute unit in sea-ice file
343  INTEGER           :: ndays_in           ! number of days
344!--- misc
345  INTEGER           :: i, j, k, l, ll     ! loop counters
346  REAL, ALLOCATABLE :: work(:,:)          ! used for extrapolation
347  CHARACTER(LEN=128):: title, mess        ! for messages
348  LOGICAL           :: extrp              ! flag for extrapolation
349  REAL              :: chmin, chmax
350  INTEGER ierr, idx
351  integer n_extrap ! number of extrapolated points
352  logical skip
353
354!------------------------------------------------------------------------------
355!---Variables depending on keyword 'mode' -------------------------------------
356  NULLIFY(champo)
357
358  SELECT CASE(mode)
359  CASE('RUG'); title='Rugosite'
360  CASE('SIC'); title='Sea-ice'
361  CASE('SST'); title='SST'
362  CASE('ALB'); title='Albedo'
363  END SELECT
364  extrp=.FALSE.; IF(PRESENT(flag).AND.mode=='SST') extrp=flag
365  idx=INDEX(fnam,'.nc')-1
366
367!--- GETTING SOME DIMENSIONAL VARIABLES FROM FILE -----------------------------
368  CALL msg(5,' Now reading file : '//TRIM(fnam))
369  CALL ncerr(NF90_OPEN(fnam, NF90_NOWRITE, ncid),fnam)
370  CALL ncerr(NF90_INQ_VARID(ncid, trim(varname), varid),fnam)
371  CALL ncerr(NF90_INQUIRE_VARIABLE(ncid, varid, dimids=dids),fnam)
372
373!--- Read unit for sea-ice variable only
374  IF (mode=='SIC') THEN
375    IF(NF90_GET_ATT(ncid, varid, 'units', unit_sic)/=NF90_NOERR) THEN
376      CALL msg(5,'No unit in sea-ice file. Take percentage as default value')
377      unit_sic='X'
378    ELSE
379      CALL msg(5,'Sea-ice cover has unit='//TRIM(unit_sic))
380    END IF
381  END IF
382
383!--- Longitude
384  CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(1), name=dnam, len=imdep),fnam)
385  ALLOCATE(dlon_ini(imdep), dlon(imdep))
386  CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam)
387  CALL ncerr(NF90_GET_VAR(ncid, varid, dlon_ini), fnam)
388  CALL msg(5,'variable '//TRIM(dnam)//' dimension ', imdep)
389
390!--- Latitude
391  CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(2), name=dnam, len=jmdep),fnam)
392  ALLOCATE(dlat_ini(jmdep), dlat(jmdep))
393  CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam)
394  CALL ncerr(NF90_GET_VAR(ncid, varid, dlat_ini), fnam)
395  CALL msg(5,'variable '//TRIM(dnam)//' dimension ', jmdep)
396
397!--- Time (variable is not needed - it is rebuilt - but calendar is)
398  CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(3), name=dnam, len=lmdep), fnam)
399  ALLOCATE(timeyear(lmdep+2))
400  CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam)
401  cal_in=' '
402  IF(NF90_GET_ATT(ncid, varid, 'calendar', cal_in)/=NF90_NOERR) THEN
403    SELECT CASE(mode)
404      CASE('RUG', 'ALB'); cal_in='360d'
405      CASE('SIC', 'SST'); cal_in='gregorian'
406    END SELECT
407  CALL msg(5,'WARNING: missing "calendar" attribute for "time" in '&
408     &//TRIM(fnam)//'. Choosing default value.')
409  END IF
410  CALL msg(5,'var, calendar, dim: '//TRIM(dnam)//' '//TRIM(cal_in), lmdep)
411 
412!--- Determining input file number of days, depending on calendar
413  ndays_in=year_len(anneeref, cal_in)
414
415!--- Rebuilding input time vector (field from input file might be unreliable)
416  IF(lmdep==12) THEN
417    timeyear=mid_month(anneeref, cal_in, ndays_in)
418    CALL msg(0,'Monthly input file(s) for '//TRIM(title)//'.')
419  ELSE IF(lmdep==ndays_in) THEN
420    timeyear=[(REAL(k)-0.5,k=0,ndays_in+1)]
421    CALL msg(0,'Daily input file (no time interpolation).')
422  ELSE
423    WRITE(mess,'(a,i3,a,i3,a)')'Mismatching input file: found',lmdep,        &
424      ' records, 12/',ndays_in,' (monthly/daily needed).'
425    CALL abort_physic('mid_months',TRIM(mess),1)
426  END IF
427
428!--- GETTING THE FIELD AND INTERPOLATING IT ----------------------------------
429  ALLOCATE(champ(imdep, jmdep), champtime(iim, jjp1, lmdep+2))
430  IF(extrp) ALLOCATE(work(imdep, jmdep))
431  CALL msg(5,'')
432  CALL msg(5,'LECTURE ET INTERPOLATION HORIZ. DE ', lmdep, ' CHAMPS.')
433  CALL ncerr(NF90_INQ_VARID(ncid, varname, varid), fnam)
434  DO l=1, lmdep
435    CALL ncerr(NF90_GET_VAR(ncid,varid,champ,[1,1,l],[imdep,jmdep,1]),fnam)
436    CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.)
437    IF(extrp) CALL extrapol(champ,imdep,jmdep,999999.,.TRUE.,.TRUE.,2,work)
438    IF(l==1) THEN
439      CALL msg(5,"----------------------------------------------------------")
440      CALL msg(5,"$$$ Interpolation barycentrique pour "//TRIM(title)//" $$$")
441      CALL msg(5,"----------------------------------------------------------")
442    END IF
443    IF(mode=='RUG') champ=LOG(champ)
444    CALL inter_barxy(dlon,dlat(:jmdep-1),champ,rlonu(:iim),rlatv,champint)
445    IF(mode=='RUG') THEN
446      champint=EXP(champint)
447      WHERE(NINT(mask)/=1) champint=0.001
448    END IF
449    champtime(:, :, l+1)=champint
450  END DO
451  CALL ncerr(NF90_CLOSE(ncid), fnam)
452
453!--- FIRST RECORD: LAST ONE OF PREVIOUS YEAR (CURRENT YEAR IF UNAVAILABLE)
454  fnam_m=fnam(1:idx)//'_m.nc'
455  IF(NF90_OPEN(fnam_m,NF90_NOWRITE,ncid)==NF90_NOERR) THEN
456    CALL msg(0,'Reading previous year file ("'//TRIM(fnam_m)//'") last record for '//TRIM(title))
457    CALL ncerr(NF90_INQ_VARID(ncid, varname, varid),fnam_m)
458    CALL ncerr(NF90_INQUIRE_VARIABLE(ncid, varid, dimids=dids),fnam_m)
459    CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(3), len=l), fnam_m)
460    CALL ncerr(NF90_GET_VAR(ncid,varid,champ,[1,1,l],[imdep,jmdep,1]),fnam_m)
461    CALL ncerr(NF90_CLOSE(ncid), fnam_m)
462    CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.)
463    IF(extrp) CALL extrapol(champ,imdep,jmdep,999999.,.TRUE.,.TRUE.,2,work)
464    IF(mode=='RUG') champ=LOG(champ)
465    CALL inter_barxy(dlon,dlat(:jmdep-1),champ,rlonu(:iim),rlatv,champint)
466    IF(mode=='RUG') THEN
467      champint=EXP(champint)
468      WHERE(NINT(mask)/=1) champint=0.001
469    END IF
470    champtime(:, :, 1)=champint
471  ELSE
472    CALL msg(0,'Using current year file ("'//TRIM(fnam)//'") last record for '//TRIM(title))
473    champtime(:, :, 1)=champtime(:, :, lmdep+1)
474  END IF
475
476!--- LAST RECORD: FIRST ONE OF NEXT YEAR (CURRENT YEAR IF UNAVAILABLE)
477  fnam_p=fnam(1:idx)//'_p.nc'
478  IF(NF90_OPEN(fnam_p,NF90_NOWRITE,ncid)==NF90_NOERR) THEN
479    CALL msg(0,'Reading previous year file ("'//TRIM(fnam_p)//'") first record for '//TRIM(title))
480    CALL ncerr(NF90_INQ_VARID(ncid, varname, varid),fnam_p)
481    CALL ncerr(NF90_GET_VAR(ncid,varid,champ,[1,1,1],[imdep,jmdep,1]),fnam_p)
482    CALL ncerr(NF90_CLOSE(ncid), fnam_p)
483    CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.)
484    IF(extrp) CALL extrapol(champ,imdep,jmdep,999999.,.TRUE.,.TRUE.,2,work)
485    IF(mode=='RUG') champ=LOG(champ)
486    CALL inter_barxy(dlon,dlat(:jmdep-1),champ,rlonu(:iim),rlatv,champint)
487    IF(mode=='RUG') THEN
488      champint=EXP(champint)
489      WHERE(NINT(mask)/=1) champint=0.001
490    END IF
491    champtime(:, :, lmdep+2)=champint
492  ELSE
493    CALL msg(0,'Using current year file ("'//TRIM(fnam)//'") first record for '//TRIM(title))
494    champtime(:, :, lmdep+2)=champtime(:, :, 2)
495  END IF
496  DEALLOCATE(dlon_ini, dlat_ini, dlon, dlat, champ)
497  IF(extrp) DEALLOCATE(work)
498
499!--- TIME INTERPOLATION ------------------------------------------------------
500  IF(prt_level>0) THEN
501     IF(ndays/=ndays_in) THEN
502        WRITE(lunout, *)
503        WRITE(lunout,*)'DIFFERENTES LONGEURS D ANNEES:'
504        WRITE(lunout,*)' Dans le fichier d entree: ',ndays_in
505        WRITE(lunout,*)' Dans le fichier de sortie: ',ndays
506     END IF
507     IF(lmdep==ndays_in) THEN
508        WRITE(lunout, *)
509        WRITE(lunout, *)'PAS D INTERPOLATION TEMPORELLE.'
510        WRITE(lunout, *)' Fichier journalier en entree.'
511     ELSE
512        WRITE(lunout, *)'INTERPOLATION TEMPORELLE.'
513        WRITE(lunout, *)' Vecteur temps en entree: ', timeyear
514        WRITE(lunout, *)' Vecteur temps en sortie de 0 a ', ndays-1
515     END IF
516  END IF
517  ALLOCATE(yder(lmdep+2), champan(iip1, jjp1, ndays))
518  IF(lmdep==ndays_in) THEN
519     champan(1:iim,:,:)=champtime
520  ELSE
521     skip = .false.
522     n_extrap = 0
523     DO j=1, jjp1
524       DO i=1, iim
525         yder = pchsp_95(timeyear, champtime(i, j, :), ibeg=2, iend=2, &
526              vc_beg=0., vc_end=0.)
527         CALL pchfe_95(timeyear, champtime(i, j, :), yder, skip, &
528              arth(0., real(ndays_in) / ndays, ndays), champan(i, j, :), ierr)
529         if (ierr < 0) stop 1
530         n_extrap = n_extrap + ierr
531       END DO
532     END DO
533     IF(n_extrap /= 0) WRITE(lunout,*) "get_2Dfield pchfe_95: n_extrap = ", n_extrap
534  END IF
535  champan(iip1, :, :)=champan(1, :, :)
536  DEALLOCATE(yder, champtime, timeyear)
537
538!--- Checking the result
539  DO j=1, jjp1
540    CALL minmax(iip1, champan(1, j, 10), chmin, chmax)
541    IF (prt_level>5) WRITE(lunout, *)' ',TRIM(title),' au temps 10 ', chmin, chmax, j
542  END DO
543
544!--- SPECIAL FILTER FOR SST: SST>271.38 --------------------------------------
545  IF(mode=='SST') THEN
546    CALL msg(5,'Filtrage de la SST: SST >= 271.38')
547    WHERE(champan<271.38) champan=271.38
548  END IF
549
550!--- SPECIAL FILTER FOR SIC: 0.0<SIC<1.0 -------------------------------------
551  IF(mode=='SIC') THEN
552    CALL msg(5,'Filtrage de la SIC: 0.0 < Sea-ice < 1.0')
553
554    IF (unit_sic=='1') THEN
555       ! Nothing to be done for sea-ice field is already in fraction of 1
556       ! This is the case for sea-ice in file cpl_atm_sic.nc
557       CALL msg(5,'Sea-ice field already in fraction of 1')
558    ELSE
559       ! Convert sea ice from percentage to fraction of 1
560       CALL msg(5,'Transformt sea-ice field from percentage to fraction of 1.')
561       champan(:, :, :)=champan(:, :, :)/100.
562    END IF
563
564    champan(iip1, :, :)=champan(1, :, :)
565    WHERE(champan>1.0) champan=1.0
566    WHERE(champan<0.0) champan=0.0
567 END IF
568
569!--- DYNAMICAL TO PHYSICAL GRID ----------------------------------------------
570  ALLOCATE(champo(klon, ndays))
571  DO k=1, ndays
572    CALL gr_dyn_fi(1, iip1, jjp1, klon, champan(1, 1, k), champo(1, k))
573  END DO
574  DEALLOCATE(champan)
575
576END SUBROUTINE get_2Dfield
577!
578!-------------------------------------------------------------------------------
579
580
581!-------------------------------------------------------------------------------
582!
583SUBROUTINE start_init_orog0(lon_in,lat_in,phis,masque)
584!
585!-------------------------------------------------------------------------------
586  USE grid_noro_m, ONLY: grid_noro0
587  IMPLICIT NONE
588!===============================================================================
589! Purpose:  Compute "phis" just like it would be in start_init_orog.
590!===============================================================================
591! Arguments:
592  REAL,             INTENT(IN)    :: lon_in(:), lat_in(:)   ! dim (iml) (jml)
593  REAL,             INTENT(INOUT) :: phis(:,:), masque(:,:) ! dim (iml,jml)
594!-------------------------------------------------------------------------------
595! Local variables:
596  CHARACTER(LEN=256) :: modname="start_init_orog0"
597  INTEGER            :: fid, llm_tmp,ttm_tmp, iml,jml, iml_rel,jml_rel, itau(1)
598  REAL               :: lev(1), date, dt, deg2rad
599  REAL, ALLOCATABLE  :: lon_rad(:), lon_ini(:), lon_rel(:,:), relief_hi(:,:)
600  REAL, ALLOCATABLE  :: lat_rad(:), lat_ini(:), lat_rel(:,:)
601!-------------------------------------------------------------------------------
602  iml=assert_eq(SIZE(lon_in),SIZE(phis,1),SIZE(masque,1),TRIM(modname)//" iml")
603  jml=assert_eq(SIZE(lat_in),SIZE(phis,2),SIZE(masque,2),TRIM(modname)//" jml")
604  IF(iml/=iip1) CALL abort_gcm(TRIM(modname),'iml/=iip1',1)
605  IF(jml/=jjp1) CALL abort_gcm(TRIM(modname),'jml/=jjp1',1)
606  pi=2.0*ASIN(1.0); deg2rad=pi/180.0
607  IF(ANY(phis/=-99999.)) RETURN                  !--- phis ALREADY KNOWN
608
609!--- HIGH RESOLUTION OROGRAPHY
610  CALL flininfo(frelf, iml_rel, jml_rel, llm_tmp, ttm_tmp, fid)
611
612  ALLOCATE(lat_rel(iml_rel,jml_rel),lon_rel(iml_rel,jml_rel))
613  CALL flinopen(frelf, .FALSE., iml_rel, jml_rel, llm_tmp, lon_rel, lat_rel,   &
614                lev, ttm_tmp, itau, date, dt, fid)
615  ALLOCATE(relief_hi(iml_rel,jml_rel))
616  CALL flinget(fid, vrel, iml_rel, jml_rel, llm_tmp, ttm_tmp, 1, 1, relief_hi)
617  CALL flinclo(fid)
618
619!--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS
620  ALLOCATE(lon_ini(iml_rel),lat_ini(jml_rel))
621  lon_ini(:)=lon_rel(:,1); IF(MAXVAL(lon_rel)>pi) lon_ini=lon_ini*deg2rad
622  lat_ini(:)=lat_rel(1,:); IF(MAXVAL(lat_rel)>pi) lat_ini=lat_ini*deg2rad
623
624!--- FIELDS ARE PROCESSED TO BE ON STANDARD ANGULAR DOMAINS
625  ALLOCATE(lon_rad(iml_rel),lat_rad(jml_rel))
626  CALL conf_dat2d(vrel, lon_ini, lat_ini, lon_rad, lat_rad, relief_hi, .FALSE.)
627  DEALLOCATE(lon_ini,lat_ini)
628
629!--- COMPUTING SURFACE GEOPOTENTIAL USING ROUTINE grid_noro0
630  WRITE(lunout,*)
631  WRITE(lunout,*)'*** Compute surface geopotential ***'
632
633!--- CALL OROGRAPHY MODULE (REDUCED VERSION) TO COMPUTE FIELDS
634  CALL grid_noro0(lon_rad, lat_rad, relief_hi, lon_in, lat_in, phis, masque)
635  phis = phis * 9.81
636  phis(iml,:) = phis(1,:)
637  DEALLOCATE(relief_hi,lon_rad,lat_rad)
638
639END SUBROUTINE start_init_orog0
640!
641!-------------------------------------------------------------------------------
642
643
644!-------------------------------------------------------------------------------
645!
646FUNCTION year_len(y,cal_in)
647!
648!-------------------------------------------------------------------------------
649  IMPLICIT NONE
650!-------------------------------------------------------------------------------
651! Arguments:
652  INTEGER                       :: year_len
653  INTEGER,           INTENT(IN) :: y
654  CHARACTER(LEN=*),  INTENT(IN) :: cal_in
655!-------------------------------------------------------------------------------
656! Local variables:
657  CHARACTER(LEN=20)             :: cal_out              ! calendar (for outputs)
658!-------------------------------------------------------------------------------
659!--- Getting the input calendar to reset at the end of the function
660  CALL ioget_calendar(cal_out)
661
662!--- Unlocking calendar and setting it to wanted one
663  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_in))
664
665!--- Getting the number of days in this year
666  year_len=ioget_year_len(y)
667
668!--- Back to original calendar
669  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_out))
670
671END FUNCTION year_len
672!
673!-------------------------------------------------------------------------------
674
675
676!-------------------------------------------------------------------------------
677!
678FUNCTION mid_month(y,cal_in,ndays_out)
679!
680!-------------------------------------------------------------------------------
681  IMPLICIT NONE
682!-------------------------------------------------------------------------------
683! Arguments:
684  INTEGER,          INTENT(IN) :: y           ! year
685  CHARACTER(LEN=*), INTENT(IN) :: cal_in      ! calendar
686  INTEGER,          INTENT(IN) :: ndays_out   ! days   number
687  REAL,          DIMENSION(14) :: mid_month   ! mid-bins times
688!-------------------------------------------------------------------------------
689! Local variables:
690  CHARACTER(LEN=99)      :: mess              ! error message
691  CHARACTER(LEN=20)      :: cal_out           ! output calendar
692  INTEGER, DIMENSION(14) :: tlen              ! months lengths (days)
693  INTEGER                :: m                 ! months counter
694  INTEGER                :: nd                ! number of days
695!-------------------------------------------------------------------------------
696  !--- Save the input calendar to reset it at the end of the function
697  CALL ioget_calendar(cal_out)
698
699  !--- Unlock calendar and set it to "cal_in"
700  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_in))
701
702  !--- Get the length of each month
703  tlen(1 )=ioget_mon_len(y-1,12)
704  DO m=1,12; tlen(m+1)=ioget_mon_len(y,m); END DO
705  tlen(14)=ioget_mon_len(y+1, 1)
706
707  !--- Mid-bins times
708  mid_month(1)=-0.5*REAL(tlen(1))
709  DO m=2,14; mid_month(m)=mid_month(m-1)+0.5*REAL(tlen(m-1)+tlen(m)); END DO
710
711  !--- Back to original calendar
712  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_out))
713
714END FUNCTION mid_month
715!
716!-------------------------------------------------------------------------------
717
718
719
720!-------------------------------------------------------------------------------
721!
722SUBROUTINE msg(lev,str1,i,str2)
723!
724!-------------------------------------------------------------------------------
725! Arguments:
726  INTEGER,                    INTENT(IN) :: lev
727  CHARACTER(LEN=*),           INTENT(IN) :: str1
728  INTEGER,          OPTIONAL, INTENT(IN) :: i
729  CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: str2
730!-------------------------------------------------------------------------------
731  IF(prt_level>lev) THEN
732    IF(PRESENT(str2)) THEN
733      WRITE(lunout,*) TRIM(str1), i, TRIM(str2)
734    ELSE IF(PRESENT(i)) THEN
735      WRITE(lunout,*) TRIM(str1), i
736    ELSE
737      WRITE(lunout,*) TRIM(str1)
738    END IF
739  END IF
740
741END SUBROUTINE msg
742!
743!-------------------------------------------------------------------------------
744
745
746!-------------------------------------------------------------------------------
747!
748SUBROUTINE ncerr(ncres,fnam)
749!
750!-------------------------------------------------------------------------------
751! Purpose: NetCDF errors handling.
752!-------------------------------------------------------------------------------
753  USE netcdf, ONLY : NF90_NOERR, NF90_STRERROR
754  IMPLICIT NONE
755!-------------------------------------------------------------------------------
756! Arguments:
757  INTEGER,          INTENT(IN) :: ncres
758  CHARACTER(LEN=*), INTENT(IN) :: fnam
759!-------------------------------------------------------------------------------
760  IF(ncres/=NF90_NOERR) THEN
761    WRITE(lunout,*)'Problem with file '//TRIM(fnam)//' in routine limit_netcdf.'
762    CALL abort_physic('limit_netcdf',NF90_STRERROR(ncres),1)
763  END IF
764
765END SUBROUTINE ncerr
766!
767!-------------------------------------------------------------------------------
768
769#endif
770! of #ifndef CPP_1D
771END SUBROUTINE limit_netcdf
772
773END MODULE limit
774!
775!*******************************************************************************
776
Note: See TracBrowser for help on using the repository browser.