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

Last change on this file since 2187 was 2187, checked in by Laurent Fairhead, 9 years ago

Merged trunk changes -r2158:2186 into testing branch.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 27.5 KB
Line 
1!
2! $Id: limit_netcdf.F90 1508 2011-04-15 13:05:34Z jghattas $
3!-------------------------------------------------------------------------------
4!
5SUBROUTINE limit_netcdf(interbar, extrap, oldice, masque)
6!
7!-------------------------------------------------------------------------------
8! Author : L. Fairhead, 27/01/94
9!-------------------------------------------------------------------------------
10! Purpose: Boundary conditions files building for new model using climatologies.
11!          Both grids have to be regular.
12!-------------------------------------------------------------------------------
13! Note: This routine is designed to work for Earth
14!-------------------------------------------------------------------------------
15! Modification history:
16!  * 23/03/1994: Z. X. Li
17!  *    09/1999: L. Fairhead (netcdf reading in LMDZ.3.3)
18!  *    07/2001: P. Le Van
19!  *    11/2009: L. Guez     (ozone day & night climatos, see etat0_netcdf.F90)
20!  *    12/2009: D. Cugnet   (f77->f90, calendars, files from coupled runs)
21!-------------------------------------------------------------------------------
22  USE control_mod
23  USE indice_sol_mod
24#ifdef CPP_EARTH
25  USE dimphy
26  USE ioipsl,             ONLY : ioget_year_len
27  USE phys_state_var_mod, ONLY : pctsrf, rlon, rlat
28  USE netcdf,             ONLY : NF90_OPEN,    NF90_CREATE,  NF90_CLOSE,       &
29                   NF90_DEF_DIM, NF90_DEF_VAR, NF90_PUT_VAR, NF90_PUT_ATT,     &
30                   NF90_NOERR,   NF90_NOWRITE, NF90_DOUBLE,  NF90_GLOBAL,      &
31                   NF90_CLOBBER, NF90_ENDDEF,  NF90_UNLIMITED, NF90_FLOAT
32  USE inter_barxy_m, only: inter_barxy
33  use netcdf95, only: nf95_def_var, nf95_put_att, nf95_put_var
34#endif
35  IMPLICIT NONE
36!-------------------------------------------------------------------------------
37! Arguments:
38#include "dimensions.h"
39#include "paramet.h"
40#include "iniprint.h"
41  LOGICAL,                    INTENT(IN) :: interbar ! barycentric interpolation
42  LOGICAL,                    INTENT(IN) :: extrap   ! SST extrapolation flag
43  LOGICAL,                    INTENT(IN) :: oldice   ! old way ice computation
44  REAL, DIMENSION(iip1,jjp1), INTENT(IN) :: masque   ! land mask
45#ifndef CPP_EARTH
46  CALL abort_gcm('limit_netcdf','Earth-specific routine, needs Earth physics',1)
47#else
48!-------------------------------------------------------------------------------
49! Local variables:
50#include "logic.h"
51#include "comvert.h"
52#include "comgeom2.h"
53#include "comconst.h"
54
55!--- INPUT NETCDF FILES NAMES --------------------------------------------------
56  CHARACTER(LEN=25) :: icefile, sstfile, dumstr
57  CHARACTER(LEN=25), PARAMETER :: famipsst ='amipbc_sst_1x1.nc        ',        &
58                                  famipsic ='amipbc_sic_1x1.nc        ',        &
59                                  fcpldsst ='cpl_atm_sst.nc           ',        &
60                                  fcpldsic ='cpl_atm_sic.nc           ',        &
61                                  fhistsst ='histmth_sst.nc           ',        &
62                                  fhistsic ='histmth_sic.nc           ',        &
63                                  frugo    ='Rugos.nc                 ',        &
64                                  falbe    ='Albedo.nc                ',        &
65                                  feraisstk='sstk.nc                  ',        &
66                                  feraici  ='ci.nc                    '
67  CHARACTER(LEN=10) :: varname
68!--- OUTPUT VARIABLES FOR NETCDF FILE ------------------------------------------
69  REAL,   DIMENSION(klon)                :: fi_ice, verif
70  REAL,   DIMENSION(:,:),   POINTER      :: phy_rug=>NULL(), phy_ice=>NULL()
71  REAL,   DIMENSION(:,:),   POINTER      :: phy_sst=>NULL(), phy_alb=>NULL()
72  REAL,   DIMENSION(:,:),   ALLOCATABLE  :: phy_bil
73  REAL,   DIMENSION(:,:,:), ALLOCATABLE  :: pctsrf_t
74  INTEGER                                :: nbad
75
76!--- VARIABLES FOR OUTPUT FILE WRITING -----------------------------------------
77  INTEGER :: ierr, nid, ndim, ntim, k
78  INTEGER, DIMENSION(2) :: dims
79  INTEGER :: id_tim,  id_SST,  id_BILS, id_RUG, id_ALB
80  INTEGER :: id_FOCE, id_FSIC, id_FTER, id_FLIC, varid_longitude, varid_latitude
81  INTEGER :: NF90_FORMAT
82  INTEGER :: ndays                   !--- Depending on the output calendar
83
84!--- INITIALIZATIONS -----------------------------------------------------------
85#ifdef NC_DOUBLE
86  NF90_FORMAT=NF90_DOUBLE
87#else
88  NF90_FORMAT=NF90_FLOAT
89#endif
90
91  pi    = 4.*ATAN(1.)
92  rad   = 6371229.
93  daysec= 86400.
94  omeg  = 2.*pi/daysec
95  g     = 9.8
96  kappa = 0.2857143
97  cpp   = 1004.70885
98  dtvr  = daysec/REAL(day_step)
99  CALL inigeom
100
101!--- Beware: anneeref (from gcm.def) is used to determine output time sampling
102  ndays=ioget_year_len(anneeref)
103
104!--- RUGOSITY TREATMENT --------------------------------------------------------
105  IF (prt_level>1) WRITE(lunout,*) 'Traitement de la rugosite'
106  varname='RUGOS'
107  CALL get_2Dfield(frugo,varname,'RUG',interbar,ndays,phy_rug,mask=masque(1:iim,:))
108
109!--- OCEAN TREATMENT -----------------------------------------------------------
110  IF (prt_level>1) WRITE(lunout,*) 'Traitement de la glace oceanique'
111
112! Input SIC file selection
113! Open file only to test if available
114  IF ( NF90_OPEN(TRIM(famipsic),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
115     icefile=TRIM(famipsic)
116     varname='sicbcs'
117  ELSE IF( NF90_OPEN(TRIM(fcpldsic),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
118     icefile=TRIM(fcpldsic)
119     varname='SIICECOV'
120  ELSE IF ( NF90_OPEN(TRIM(fhistsic),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
121     icefile=TRIM(fhistsic)
122     varname='pourc_sic'
123  ELSE IF ( NF90_OPEN(TRIM(feraici),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
124     icefile=TRIM(feraici)
125     varname='ci'
126  ELSE
127     WRITE(lunout,*) 'ERROR! No sea-ice input file was found.'
128
129     WRITE(lunout,*) 'One of following files must be availible : ',trim(famipsic),', ',trim(fcpldsic),', ', &
130                      trim(fhistsic), trim(feraici)
131
132     CALL abort_gcm('limit_netcdf','No sea-ice file was found',1)
133  END IF
134  ierr=NF90_CLOSE(nid)
135  IF (prt_level>=0) WRITE(lunout,*)'Pour la glace de mer a ete choisi le fichier '//TRIM(icefile)
136
137  CALL get_2Dfield(icefile,varname, 'SIC',interbar,ndays,phy_ice,flag=oldice)
138
139  ALLOCATE(pctsrf_t(klon,nbsrf,ndays))
140  DO k=1,ndays
141     fi_ice=phy_ice(:,k)
142     WHERE(fi_ice>=1.0  ) fi_ice=1.0
143     WHERE(fi_ice<EPSFRA) fi_ice=0.0
144     pctsrf_t(:,is_ter,k)=pctsrf(:,is_ter)       ! land soil
145     pctsrf_t(:,is_lic,k)=pctsrf(:,is_lic)       ! land ice
146     IF (icefile==trim(fcpldsic)) THEN           ! SIC=pICE*(1-LIC-TER)
147        pctsrf_t(:,is_sic,k)=fi_ice(:)*(1.-pctsrf(:,is_lic)-pctsrf(:,is_ter))
148     ELSE IF (icefile==trim(fhistsic)) THEN      ! SIC=pICE
149        pctsrf_t(:,is_sic,k)=fi_ice(:)
150     ELSE ! icefile==famipsic                    ! SIC=pICE-LIC
151        pctsrf_t(:,is_sic,k)=fi_ice-pctsrf_t(:,is_lic,k)
152     END IF
153     WHERE(pctsrf_t(:,is_sic,k)<=0) pctsrf_t(:,is_sic,k)=0.
154     WHERE(1.0-zmasq<EPSFRA)
155        pctsrf_t(:,is_sic,k)=0.0
156        pctsrf_t(:,is_oce,k)=0.0
157     ELSEWHERE
158        WHERE(pctsrf_t(:,is_sic,k)>=1.0-zmasq)
159           pctsrf_t(:,is_sic,k)=1.0-zmasq
160           pctsrf_t(:,is_oce,k)=0.0
161        ELSEWHERE
162           pctsrf_t(:,is_oce,k)=1.0-zmasq-pctsrf_t(:,is_sic,k)
163           WHERE(pctsrf_t(:,is_oce,k)<EPSFRA)
164              pctsrf_t(:,is_oce,k)=0.0
165              pctsrf_t(:,is_sic,k)=1.0-zmasq
166           END WHERE
167        END WHERE
168     END WHERE
169     nbad=COUNT(pctsrf_t(:,is_oce,k)<0.0)
170     IF(nbad>0) WRITE(lunout,*) 'pb sous maille pour nb point = ',nbad
171     nbad=COUNT(abs(sum(pctsrf_t(:,:,k),dim=2)-1.0)>EPSFRA)
172     IF(nbad>0) WRITE(lunout,*) 'pb sous surface pour nb points = ',nbad
173  END DO
174  DEALLOCATE(phy_ice)
175
176!--- SST TREATMENT -------------------------------------------------------------
177  IF (prt_level>1) WRITE(lunout,*) 'Traitement de la sst'
178
179! Input SST file selection
180! Open file only to test if available
181  IF ( NF90_OPEN(TRIM(famipsst),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
182     sstfile=TRIM(famipsst)
183     varname='tosbcs'
184  ELSE IF ( NF90_OPEN(TRIM(fcpldsst),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
185     sstfile=TRIM(fcpldsst)
186     varname='SISUTESW'
187  ELSE IF ( NF90_OPEN(TRIM(fhistsst),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
188     sstfile=TRIM(fhistsst)
189     varname='tsol_oce'
190  ELSE IF ( NF90_OPEN(TRIM(feraisstk),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
191     sstfile=TRIM(feraisstk)
192     varname='sstk'
193  ELSE
194     WRITE(lunout,*) 'ERROR! No sst input file was found.'
195     WRITE(lunout,*) 'One of following files must be availible : ',trim(famipsst),trim(fcpldsst),trim(fhistsst),trim(feraisstk)
196     CALL abort_gcm('limit_netcdf','No sst file was found',1)
197  END IF
198  ierr=NF90_CLOSE(nid)
199  IF (prt_level>=0) WRITE(lunout,*)'Pour la temperature de mer a ete choisi un fichier '//TRIM(sstfile)
200
201  CALL get_2Dfield(sstfile,varname,'SST',interbar,ndays,phy_sst,flag=extrap)
202
203!--- ALBEDO TREATMENT ----------------------------------------------------------
204  IF (prt_level>1) WRITE(lunout,*) 'Traitement de l albedo'
205  varname='ALBEDO'
206  CALL get_2Dfield(falbe,varname,'ALB',interbar,ndays,phy_alb)
207
208!--- REFERENCE GROUND HEAT FLUX TREATMENT --------------------------------------
209  ALLOCATE(phy_bil(klon,ndays)); phy_bil=0.0
210
211!--- OUTPUT FILE WRITING -------------------------------------------------------
212  IF (prt_level>5) WRITE(lunout,*) 'Ecriture du fichier limit : debut'
213
214  !--- File creation
215  ierr=NF90_CREATE("limit.nc",NF90_CLOBBER,nid)
216  ierr=NF90_PUT_ATT(nid,NF90_GLOBAL,"title","Fichier conditions aux limites")
217
218  !--- Dimensions creation
219  ierr=NF90_DEF_DIM(nid,"points_physiques",klon,ndim)
220  ierr=NF90_DEF_DIM(nid,"time",NF90_UNLIMITED,ntim)
221
222  dims=(/ndim,ntim/)
223
224  !--- Variables creation
225  ierr=NF90_DEF_VAR(nid,"TEMPS",NF90_FORMAT,(/ntim/),id_tim)
226  ierr=NF90_DEF_VAR(nid,"FOCE", NF90_FORMAT,dims,id_FOCE)
227  ierr=NF90_DEF_VAR(nid,"FSIC", NF90_FORMAT,dims,id_FSIC)
228  ierr=NF90_DEF_VAR(nid,"FTER", NF90_FORMAT,dims,id_FTER)
229  ierr=NF90_DEF_VAR(nid,"FLIC", NF90_FORMAT,dims,id_FLIC)
230  ierr=NF90_DEF_VAR(nid,"SST",  NF90_FORMAT,dims,id_SST)
231  ierr=NF90_DEF_VAR(nid,"BILS", NF90_FORMAT,dims,id_BILS)
232  ierr=NF90_DEF_VAR(nid,"ALB",  NF90_FORMAT,dims,id_ALB)
233  ierr=NF90_DEF_VAR(nid,"RUG",  NF90_FORMAT,dims,id_RUG)
234  call nf95_def_var(nid, "longitude", NF90_FLOAT, ndim, varid_longitude)
235  call nf95_def_var(nid, "latitude", NF90_FLOAT, ndim, varid_latitude)
236
237  !--- Attributes creation
238  ierr=NF90_PUT_ATT(nid,id_tim, "title","Jour dans l annee")
239  ierr=NF90_PUT_ATT(nid,id_FOCE,"title","Fraction ocean")
240  ierr=NF90_PUT_ATT(nid,id_FSIC,"title","Fraction glace de mer")
241  ierr=NF90_PUT_ATT(nid,id_FTER,"title","Fraction terre")
242  ierr=NF90_PUT_ATT(nid,id_FLIC,"title","Fraction land ice")
243  ierr=NF90_PUT_ATT(nid,id_SST ,"title","Temperature superficielle de la mer")
244  ierr=NF90_PUT_ATT(nid,id_BILS,"title","Reference flux de chaleur au sol")
245  ierr=NF90_PUT_ATT(nid,id_ALB, "title","Albedo a la surface")
246  ierr=NF90_PUT_ATT(nid,id_RUG, "title","Rugosite")
247
248  call nf95_put_att(nid, varid_longitude, "standard_name", "longitude")
249  call nf95_put_att(nid, varid_longitude, "units", "degrees_east")
250
251  call nf95_put_att(nid, varid_latitude, "standard_name", "latitude")
252  call nf95_put_att(nid, varid_latitude, "units", "degrees_north")
253
254  ierr=NF90_ENDDEF(nid)
255
256  !--- Variables saving
257  ierr=NF90_PUT_VAR(nid,id_tim,(/(REAL(k),k=1,ndays)/))
258  ierr=NF90_PUT_VAR(nid,id_FOCE,pctsrf_t(:,is_oce,:),(/1,1/),(/klon,ndays/))
259  ierr=NF90_PUT_VAR(nid,id_FSIC,pctsrf_t(:,is_sic,:),(/1,1/),(/klon,ndays/))
260  ierr=NF90_PUT_VAR(nid,id_FTER,pctsrf_t(:,is_ter,:),(/1,1/),(/klon,ndays/))
261  ierr=NF90_PUT_VAR(nid,id_FLIC,pctsrf_t(:,is_lic,:),(/1,1/),(/klon,ndays/))
262  ierr=NF90_PUT_VAR(nid,id_SST ,phy_sst(:,:),(/1,1/),(/klon,ndays/))
263  ierr=NF90_PUT_VAR(nid,id_BILS,phy_bil(:,:),(/1,1/),(/klon,ndays/))
264  ierr=NF90_PUT_VAR(nid,id_ALB ,phy_alb(:,:),(/1,1/),(/klon,ndays/))
265  ierr=NF90_PUT_VAR(nid,id_RUG ,phy_rug(:,:),(/1,1/),(/klon,ndays/))
266  call nf95_put_var(nid, varid_longitude, rlon)
267  call nf95_put_var(nid, varid_latitude, rlat)
268
269  ierr=NF90_CLOSE(nid)
270
271  IF (prt_level>5) WRITE(lunout,*) 'Ecriture du fichier limit : fin'
272
273  DEALLOCATE(pctsrf_t,phy_sst,phy_bil,phy_alb,phy_rug)
274
275
276!===============================================================================
277!
278  CONTAINS
279!
280!===============================================================================
281
282
283!-------------------------------------------------------------------------------
284!
285SUBROUTINE get_2Dfield(fnam, varname, mode, ibar, ndays, champo, flag, mask)
286!
287!-----------------------------------------------------------------------------
288! Comments:
289!   There are two assumptions concerning the NetCDF files, that are satisfied
290!   with files that are conforming NC convention:
291!     1) The last dimension of the variables used is the time record.
292!     2) Dimensional variables have the same names as corresponding dimensions.
293!-----------------------------------------------------------------------------
294  USE netcdf, ONLY: NF90_OPEN, NF90_INQ_VARID, NF90_INQUIRE_VARIABLE, &
295       NF90_CLOSE, NF90_INQ_DIMID, NF90_INQUIRE_DIMENSION, NF90_GET_VAR, &
296       NF90_GET_ATT
297  USE dimphy, ONLY : klon
298  USE phys_state_var_mod, ONLY : pctsrf
299  USE control_mod
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#include "iniprint.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  LOGICAL,           INTENT(IN)     :: ibar     ! interp on pressure levels
316  INTEGER,           INTENT(IN)     :: ndays    ! current year number of days
317  REAL,    POINTER,  DIMENSION(:, :) :: champo  ! output field = f(t)
318  LOGICAL, OPTIONAL, INTENT(IN)     :: flag     ! extrapol. (SST) old ice (SIC)
319  REAL,    OPTIONAL, DIMENSION(iim, jjp1), INTENT(IN) :: mask
320!------------------------------------------------------------------------------
321! Local variables:
322!--- NetCDF
323  INTEGER :: ncid, varid                  ! NetCDF identifiers
324  CHARACTER(LEN=30)               :: dnam       ! dimension name
325!--- dimensions
326  INTEGER,           DIMENSION(4) :: dids       ! NetCDF dimensions identifiers
327  REAL, ALLOCATABLE, DIMENSION(:) :: dlon_ini   ! initial longitudes vector
328  REAL, ALLOCATABLE, DIMENSION(:) :: dlat_ini   ! initial latitudes  vector
329  REAL, POINTER,     DIMENSION(:) :: dlon, dlat ! reordered lon/lat  vectors
330!--- fields
331  INTEGER :: imdep, jmdep, lmdep                ! dimensions of 'champ'
332  REAL, ALLOCATABLE, DIMENSION(:, :) :: champ   ! wanted field on initial grid
333  REAL, ALLOCATABLE, DIMENSION(:) :: yder, timeyear
334  REAL,              DIMENSION(iim, jjp1) :: champint   ! interpolated field
335  REAL, ALLOCATABLE, DIMENSION(:, :, :)    :: champtime
336  REAL, ALLOCATABLE, DIMENSION(:, :, :)    :: champan
337!--- input files
338  CHARACTER(LEN=20)                 :: cal_in   ! calendar
339  CHARACTER(LEN=20)                 :: unit_sic ! attribute unit in sea-ice file
340  INTEGER                           :: ndays_in ! number of days
341!--- misc
342  INTEGER :: i, j, k, l                         ! loop counters
343  REAL, ALLOCATABLE, DIMENSION(:, :) :: work     ! used for extrapolation
344  CHARACTER(LEN=25)                 :: title    ! for messages
345  LOGICAL                           :: extrp    ! flag for extrapolation
346  LOGICAL                           :: oldice   ! flag for old way ice computation
347  REAL                              :: chmin, chmax
348  INTEGER ierr
349  integer n_extrap ! number of extrapolated points
350  logical skip
351
352!------------------------------------------------------------------------------
353!---Variables depending on keyword 'mode' -------------------------------------
354  NULLIFY(champo)
355
356  SELECT CASE(mode)
357  CASE('RUG'); title='Rugosite'
358  CASE('SIC'); title='Sea-ice'
359  CASE('SST'); title='SST'
360  CASE('ALB'); title='Albedo'
361  END SELECT
362 
363
364  extrp=.FALSE.
365  oldice=.FALSE.
366  IF ( PRESENT(flag) ) THEN
367    IF ( flag .AND. mode=='SST' ) extrp=.TRUE.
368    IF ( flag .AND. mode=='SIC' ) oldice=.TRUE.
369  END IF
370
371!--- GETTING SOME DIMENSIONAL VARIABLES FROM FILE -----------------------------
372  IF (prt_level>5) WRITE(lunout,*) ' Now reading file : ',fnam
373  ierr=NF90_OPEN(fnam, NF90_NOWRITE, ncid);             CALL ncerr(ierr, fnam)
374  ierr=NF90_INQ_VARID(ncid, trim(varname), varid);            CALL ncerr(ierr, fnam)
375  ierr=NF90_INQUIRE_VARIABLE(ncid, varid, dimids=dids); CALL ncerr(ierr, fnam)
376
377!--- Read unit for sea-ice variable only
378  IF (mode=='SIC') THEN
379     ierr=NF90_GET_ATT(ncid, varid, 'units', unit_sic)
380     IF(ierr/=NF90_NOERR) THEN
381        IF (prt_level>5) WRITE(lunout,*) 'No unit was given in sea-ice file. Take percentage as default value'
382        unit_sic='X'
383     ELSE
384        IF (prt_level>5) WRITE(lunout,*) ' Sea-ice cover has unit=',unit_sic
385     END IF
386  END IF
387
388!--- Longitude
389  ierr=NF90_INQUIRE_DIMENSION(ncid, dids(1), name=dnam, len=imdep)
390  CALL ncerr(ierr, fnam); ALLOCATE(dlon_ini(imdep), dlon(imdep))
391  ierr=NF90_INQ_VARID(ncid, dnam, varid);                CALL ncerr(ierr, fnam)
392  ierr=NF90_GET_VAR(ncid, varid, dlon_ini);              CALL ncerr(ierr, fnam)
393  IF (prt_level>5) WRITE(lunout, *) 'variable ', dnam, 'dimension ', imdep
394
395!--- Latitude
396  ierr=NF90_INQUIRE_DIMENSION(ncid, dids(2), name=dnam, len=jmdep)
397  CALL ncerr(ierr, fnam); ALLOCATE(dlat_ini(jmdep), dlat(jmdep))
398  ierr=NF90_INQ_VARID(ncid, dnam, varid);                CALL ncerr(ierr, fnam)
399  ierr=NF90_GET_VAR(ncid, varid, dlat_ini);              CALL ncerr(ierr, fnam)
400  IF (prt_level>5) WRITE(lunout, *) 'variable ', dnam, 'dimension ', jmdep
401
402!--- Time (variable is not needed - it is rebuilt - but calendar is)
403  ierr=NF90_INQUIRE_DIMENSION(ncid, dids(3), name=dnam, len=lmdep)
404  CALL ncerr(ierr, fnam); ALLOCATE(timeyear(lmdep))
405  ierr=NF90_INQ_VARID(ncid, dnam, varid);                CALL ncerr(ierr, fnam)
406  cal_in=' '
407  ierr=NF90_GET_ATT(ncid, varid, 'calendar', cal_in)
408  IF(ierr/=NF90_NOERR) THEN
409    SELECT CASE(mode)
410      CASE('RUG', 'ALB'); cal_in='360d'
411      CASE('SIC', 'SST'); cal_in='gregorian'
412    END SELECT
413    IF (prt_level>5) WRITE(lunout, *)'ATTENTION: variable "time" sans attribut "calendrier" ' &
414         // 'dans '//TRIM(fnam)//'. On choisit la valeur par defaut.'
415  END IF
416  IF (prt_level>5) WRITE(lunout, *) 'variable ', dnam, 'dimension ', lmdep, 'calendrier ', &
417       cal_in
418
419 
420!--- CONSTRUCTING THE INPUT TIME VECTOR FOR INTERPOLATION --------------------
421  !--- Determining input file number of days, depending on calendar
422  ndays_in=year_len(anneeref, cal_in)
423
424!--- Time vector reconstruction (time vector from file is not trusted)
425!--- If input records are not monthly, time sampling has to be constant !
426  timeyear=mid_months(anneeref, cal_in, lmdep)
427  IF (lmdep /= 12) WRITE(lunout,*) 'Note : les fichiers de ', TRIM(mode), &
428       ' ne comportent pas 12, mais ', lmdep, ' enregistrements.'
429
430!--- GETTING THE FIELD AND INTERPOLATING IT ----------------------------------
431  ALLOCATE(champ(imdep, jmdep), champtime(iim, jjp1, lmdep))
432  IF(extrp) ALLOCATE(work(imdep, jmdep))
433
434  IF (prt_level>5) WRITE(lunout, *)
435  IF (prt_level>5) WRITE(lunout,*)'LECTURE ET INTERPOLATION HORIZ. DE ', lmdep, ' CHAMPS.'
436  ierr=NF90_INQ_VARID(ncid, varname, varid);             CALL ncerr(ierr, fnam)
437  DO l=1, lmdep
438    ierr=NF90_GET_VAR(ncid, varid, champ, (/1, 1, l/), (/imdep, jmdep, 1/))
439    CALL ncerr(ierr, fnam)
440    CALL conf_dat2d(title, imdep, jmdep, dlon_ini, dlat_ini, dlon, dlat, &
441         champ, ibar)
442
443    IF (extrp) CALL extrapol(champ, imdep, jmdep, 999999., .TRUE., .TRUE., 2, &
444         work)
445
446    IF(ibar .AND. .NOT.oldice) THEN
447      IF(l==1 .AND. prt_level>5) THEN
448        WRITE(lunout, *) '-------------------------------------------------------------------------'
449        WRITE(lunout, *) 'Utilisation de l''interpolation barycentrique pour ',TRIM(title),' $$$'
450        WRITE(lunout, *) '-------------------------------------------------------------------------'
451      END IF
452      IF(mode=='RUG') champ=LOG(champ)
453      CALL inter_barxy(dlon, dlat(:jmdep-1), champ, rlonu(:iim),     &
454                         rlatv, champint)
455      IF(mode=='RUG') THEN
456        champint=EXP(champint)
457        WHERE(NINT(mask)/=1) champint=0.001
458      END IF
459    ELSE
460      SELECT CASE(mode)
461        CASE('RUG');       CALL rugosite(imdep, jmdep, dlon, dlat, champ,    &
462                                    iim, jjp1, rlonv, rlatu, champint, mask)
463        CASE('SIC');       CALL sea_ice (imdep, jmdep, dlon, dlat, champ,    &
464                                    iim, jjp1, rlonv, rlatu, champint)
465        CASE('SST', 'ALB'); CALL grille_m(imdep, jmdep, dlon, dlat, champ,    &
466                                    iim, jjp1, rlonv, rlatu, champint)
467      END SELECT
468    END IF
469    champtime(:, :, l)=champint
470  END DO
471  ierr=NF90_CLOSE(ncid);                                 CALL ncerr(ierr, fnam)
472
473  DEALLOCATE(dlon_ini, dlat_ini, dlon, dlat, champ)
474  IF(extrp) DEALLOCATE(work)
475
476!--- TIME INTERPOLATION ------------------------------------------------------
477  IF (prt_level>5) THEN
478     WRITE(lunout, *)
479     WRITE(lunout, *)'INTERPOLATION TEMPORELLE.'
480     WRITE(lunout, *)' Vecteur temps en entree: ', timeyear
481     WRITE(lunout, *)' Vecteur temps en sortie de 0 a ', ndays
482  END IF
483
484  ALLOCATE(yder(lmdep), champan(iip1, jjp1, ndays))
485  skip = .false.
486  n_extrap = 0
487  DO j=1, jjp1
488    DO i=1, iim
489      yder = pchsp_95(timeyear, champtime(i, j, :), ibeg=2, iend=2, &
490           vc_beg=0., vc_end=0.)
491      CALL pchfe_95(timeyear, champtime(i, j, :), yder, skip, &
492           arth(0., real(ndays_in) / ndays, ndays), champan(i, j, :), ierr)
493      if (ierr < 0) stop 1
494      n_extrap = n_extrap + ierr
495    END DO
496  END DO
497  if (n_extrap /= 0) then
498     WRITE(lunout,*) "get_2Dfield pchfe_95: n_extrap = ", n_extrap
499  end if
500  champan(iip1, :, :)=champan(1, :, :)
501  DEALLOCATE(yder, champtime, timeyear)
502
503!--- Checking the result
504  DO j=1, jjp1
505    CALL minmax(iip1, champan(1, j, 10), chmin, chmax)
506    IF (prt_level>5) WRITE(lunout, *)' ',TRIM(title),' au temps 10 ', chmin, chmax, j
507  END DO
508
509!--- SPECIAL FILTER FOR SST: SST>271.38 --------------------------------------
510  IF(mode=='SST') THEN
511    IF (prt_level>5) WRITE(lunout, *) 'Filtrage de la SST: SST >= 271.38'
512    WHERE(champan<271.38) champan=271.38
513  END IF
514
515!--- SPECIAL FILTER FOR SIC: 0.0<SIC<1.0 -------------------------------------
516  IF(mode=='SIC') THEN
517    IF (prt_level>5) WRITE(lunout, *) 'Filtrage de la SIC: 0.0 < Sea-ice < 1.0'
518
519    IF (unit_sic=='1') THEN
520       ! Nothing to be done for sea-ice field is already in fraction of 1
521       ! This is the case for sea-ice in file cpl_atm_sic.nc
522       IF (prt_level>5) WRITE(lunout,*) 'Sea-ice field already in fraction of 1'
523    ELSE
524       ! Convert sea ice from percentage to fraction of 1
525       IF (prt_level>5) WRITE(lunout,*) 'Transformt sea-ice field from percentage to fraction of 1.'
526       champan(:, :, :)=champan(:, :, :)/100.
527    END IF
528
529    champan(iip1, :, :)=champan(1, :, :)
530    WHERE(champan>1.0) champan=1.0
531    WHERE(champan<0.0) champan=0.0
532 END IF
533
534!--- DYNAMICAL TO PHYSICAL GRID ----------------------------------------------
535  ALLOCATE(champo(klon, ndays))
536  DO k=1, ndays
537    CALL gr_dyn_fi(1, iip1, jjp1, klon, champan(1, 1, k), champo(1, k))
538  END DO
539  DEALLOCATE(champan)
540
541END SUBROUTINE get_2Dfield
542!
543!-------------------------------------------------------------------------------
544
545
546
547!-------------------------------------------------------------------------------
548!
549FUNCTION year_len(y,cal_in)
550!
551!-------------------------------------------------------------------------------
552  USE ioipsl, ONLY : ioget_calendar,ioconf_calendar,lock_calendar,ioget_year_len
553  IMPLICIT NONE
554!-------------------------------------------------------------------------------
555! Arguments:
556  INTEGER                       :: year_len
557  INTEGER,           INTENT(IN) :: y
558  CHARACTER(LEN=*),  INTENT(IN) :: cal_in
559!-------------------------------------------------------------------------------
560! Local variables:
561  CHARACTER(LEN=20)             :: cal_out              ! calendar (for outputs)
562!-------------------------------------------------------------------------------
563!--- Getting the input calendar to reset at the end of the function
564  CALL ioget_calendar(cal_out)
565
566!--- Unlocking calendar and setting it to wanted one
567  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_in))
568
569!--- Getting the number of days in this year
570  year_len=ioget_year_len(y)
571
572!--- Back to original calendar
573  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_out))
574
575END FUNCTION year_len
576!
577!-------------------------------------------------------------------------------
578
579
580!-------------------------------------------------------------------------------
581!
582FUNCTION mid_months(y,cal_in,nm)
583!
584!-------------------------------------------------------------------------------
585  USE ioipsl, ONLY : ioget_calendar,ioconf_calendar,lock_calendar,ioget_mon_len
586  IMPLICIT NONE
587!-------------------------------------------------------------------------------
588! Arguments:
589  INTEGER,                INTENT(IN) :: y               ! year
590  CHARACTER(LEN=*),       INTENT(IN) :: cal_in          ! calendar
591  INTEGER,                INTENT(IN) :: nm              ! months/year number
592  REAL,    DIMENSION(nm)             :: mid_months      ! mid-month times
593!-------------------------------------------------------------------------------
594! Local variables:
595  CHARACTER(LEN=99)                  :: mess            ! error message
596  CHARACTER(LEN=20)                  :: cal_out         ! calendar (for outputs)
597  INTEGER, DIMENSION(nm)             :: mnth            ! months lengths (days)
598  INTEGER                            :: m               ! months counter
599  INTEGER                            :: nd              ! number of days
600!-------------------------------------------------------------------------------
601  nd=year_len(y,cal_in)
602
603  IF(nm==12) THEN
604
605  !--- Getting the input calendar to reset at the end of the function
606    CALL ioget_calendar(cal_out)
607
608  !--- Unlocking calendar and setting it to wanted one
609    CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_in))
610
611  !--- Getting the length of each month
612    DO m=1,nm; mnth(m)=ioget_mon_len(y,m); END DO
613
614  !--- Back to original calendar
615    CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_out))
616
617  ELSE IF(MODULO(nd,nm)/=0) THEN
618    WRITE(mess,'(a,i3,a,i3,a)')'Unconsistent calendar: ',nd,' days/year, but ',&
619      nm,' months/year. Months number should divide days number.'
620    CALL abort_gcm('mid_months',TRIM(mess),1)
621
622  ELSE
623    mnth=(/(m,m=1,nm,nd/nm)/)
624  END IF
625
626!--- Mid-months times
627  mid_months(1)=0.5*REAL(mnth(1))
628  DO k=2,nm
629    mid_months(k)=mid_months(k-1)+0.5*REAL(mnth(k-1)+mnth(k))
630  END DO
631
632END FUNCTION mid_months
633!
634!-------------------------------------------------------------------------------
635
636
637!-------------------------------------------------------------------------------
638!
639SUBROUTINE ncerr(ncres,fnam)
640!
641!-------------------------------------------------------------------------------
642! Purpose: NetCDF errors handling.
643!-------------------------------------------------------------------------------
644  USE netcdf, ONLY : NF90_NOERR, NF90_STRERROR
645  IMPLICIT NONE
646!-------------------------------------------------------------------------------
647! Arguments:
648  INTEGER,          INTENT(IN) :: ncres
649  CHARACTER(LEN=*), INTENT(IN) :: fnam
650!-------------------------------------------------------------------------------
651#include "iniprint.h"
652  IF(ncres/=NF90_NOERR) THEN
653    WRITE(lunout,*)'Problem with file '//TRIM(fnam)//' in routine limit_netcdf.'
654    CALL abort_gcm('limit_netcdf',NF90_STRERROR(ncres),1)
655  END IF
656
657END SUBROUTINE ncerr
658!
659!-------------------------------------------------------------------------------
660
661#endif
662! of #ifdef CPP_EARTH
663
664END SUBROUTINE limit_netcdf
Note: See TracBrowser for help on using the repository browser.