source: LMDZ5/trunk/libf/phylmd/limit_netcdf.F90 @ 2192

Last change on this file since 2192 was 2159, checked in by jyg, 10 years ago

1/ Splitting of the boundary layer : the climbing down and up of Pbl_surface is
split between the off-wake and wake regions ; the thermal scheme is applied
only to the off-wake region.
2/ Elimination of wake_scal and calwake_scal.

  • 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     WRITE(lunout,*) 'One of following files must be availible : ',trim(famipsic),', ',trim(fcpldsic),', ', &
129                      trim(fhistsic), trim(feraici)
130     CALL abort_gcm('limit_netcdf','No sea-ice file was found',1)
131  END IF
132  ierr=NF90_CLOSE(nid)
133  IF (prt_level>=0) WRITE(lunout,*)'Pour la glace de mer a ete choisi le fichier '//TRIM(icefile)
134
135  CALL get_2Dfield(icefile,varname, 'SIC',interbar,ndays,phy_ice,flag=oldice)
136
137  ALLOCATE(pctsrf_t(klon,nbsrf,ndays))
138  DO k=1,ndays
139     fi_ice=phy_ice(:,k)
140     WHERE(fi_ice>=1.0  ) fi_ice=1.0
141     WHERE(fi_ice<EPSFRA) fi_ice=0.0
142     pctsrf_t(:,is_ter,k)=pctsrf(:,is_ter)       ! land soil
143     pctsrf_t(:,is_lic,k)=pctsrf(:,is_lic)       ! land ice
144     IF (icefile==trim(fcpldsic)) THEN           ! SIC=pICE*(1-LIC-TER)
145        pctsrf_t(:,is_sic,k)=fi_ice(:)*(1.-pctsrf(:,is_lic)-pctsrf(:,is_ter))
146     ELSE IF (icefile==trim(fhistsic)) THEN      ! SIC=pICE
147        pctsrf_t(:,is_sic,k)=fi_ice(:)
148     ELSE ! icefile==famipsic                    ! SIC=pICE-LIC
149        pctsrf_t(:,is_sic,k)=fi_ice-pctsrf_t(:,is_lic,k)
150     END IF
151     WHERE(pctsrf_t(:,is_sic,k)<=0) pctsrf_t(:,is_sic,k)=0.
152     WHERE(1.0-zmasq<EPSFRA)
153        pctsrf_t(:,is_sic,k)=0.0
154        pctsrf_t(:,is_oce,k)=0.0
155     ELSEWHERE
156        WHERE(pctsrf_t(:,is_sic,k)>=1.0-zmasq)
157           pctsrf_t(:,is_sic,k)=1.0-zmasq
158           pctsrf_t(:,is_oce,k)=0.0
159        ELSEWHERE
160           pctsrf_t(:,is_oce,k)=1.0-zmasq-pctsrf_t(:,is_sic,k)
161           WHERE(pctsrf_t(:,is_oce,k)<EPSFRA)
162              pctsrf_t(:,is_oce,k)=0.0
163              pctsrf_t(:,is_sic,k)=1.0-zmasq
164           END WHERE
165        END WHERE
166     END WHERE
167     nbad=COUNT(pctsrf_t(:,is_oce,k)<0.0)
168     IF(nbad>0) WRITE(lunout,*) 'pb sous maille pour nb point = ',nbad
169     nbad=COUNT(abs(sum(pctsrf_t(:,:,k),dim=2)-1.0)>EPSFRA)
170     IF(nbad>0) WRITE(lunout,*) 'pb sous surface pour nb points = ',nbad
171  END DO
172  DEALLOCATE(phy_ice)
173
174!--- SST TREATMENT -------------------------------------------------------------
175  IF (prt_level>1) WRITE(lunout,*) 'Traitement de la sst'
176
177! Input SST file selection
178! Open file only to test if available
179  IF ( NF90_OPEN(TRIM(famipsst),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
180     sstfile=TRIM(famipsst)
181     varname='tosbcs'
182  ELSE IF ( NF90_OPEN(TRIM(fcpldsst),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
183     sstfile=TRIM(fcpldsst)
184     varname='SISUTESW'
185  ELSE IF ( NF90_OPEN(TRIM(fhistsst),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
186     sstfile=TRIM(fhistsst)
187     varname='tsol_oce'
188  ELSE IF ( NF90_OPEN(TRIM(feraisstk),NF90_NOWRITE,nid)==NF90_NOERR ) THEN
189     sstfile=TRIM(feraisstk)
190     varname='sstk'
191  ELSE
192     WRITE(lunout,*) 'ERROR! No sst input file was found.'
193     WRITE(lunout,*) 'One of following files must be availible : ',trim(famipsst),trim(fcpldsst),trim(fhistsst),trim(feraisstk)
194     CALL abort_gcm('limit_netcdf','No sst file was found',1)
195  END IF
196  ierr=NF90_CLOSE(nid)
197  IF (prt_level>=0) WRITE(lunout,*)'Pour la temperature de mer a ete choisi un fichier '//TRIM(sstfile)
198
199  CALL get_2Dfield(sstfile,varname,'SST',interbar,ndays,phy_sst,flag=extrap)
200
201!--- ALBEDO TREATMENT ----------------------------------------------------------
202  IF (prt_level>1) WRITE(lunout,*) 'Traitement de l albedo'
203  varname='ALBEDO'
204  CALL get_2Dfield(falbe,varname,'ALB',interbar,ndays,phy_alb)
205
206!--- REFERENCE GROUND HEAT FLUX TREATMENT --------------------------------------
207  ALLOCATE(phy_bil(klon,ndays)); phy_bil=0.0
208
209!--- OUTPUT FILE WRITING -------------------------------------------------------
210  IF (prt_level>5) WRITE(lunout,*) 'Ecriture du fichier limit : debut'
211
212  !--- File creation
213  ierr=NF90_CREATE("limit.nc",NF90_CLOBBER,nid)
214  ierr=NF90_PUT_ATT(nid,NF90_GLOBAL,"title","Fichier conditions aux limites")
215
216  !--- Dimensions creation
217  ierr=NF90_DEF_DIM(nid,"points_physiques",klon,ndim)
218  ierr=NF90_DEF_DIM(nid,"time",NF90_UNLIMITED,ntim)
219
220  dims=(/ndim,ntim/)
221
222  !--- Variables creation
223  ierr=NF90_DEF_VAR(nid,"TEMPS",NF90_FORMAT,(/ntim/),id_tim)
224  ierr=NF90_DEF_VAR(nid,"FOCE", NF90_FORMAT,dims,id_FOCE)
225  ierr=NF90_DEF_VAR(nid,"FSIC", NF90_FORMAT,dims,id_FSIC)
226  ierr=NF90_DEF_VAR(nid,"FTER", NF90_FORMAT,dims,id_FTER)
227  ierr=NF90_DEF_VAR(nid,"FLIC", NF90_FORMAT,dims,id_FLIC)
228  ierr=NF90_DEF_VAR(nid,"SST",  NF90_FORMAT,dims,id_SST)
229  ierr=NF90_DEF_VAR(nid,"BILS", NF90_FORMAT,dims,id_BILS)
230  ierr=NF90_DEF_VAR(nid,"ALB",  NF90_FORMAT,dims,id_ALB)
231  ierr=NF90_DEF_VAR(nid,"RUG",  NF90_FORMAT,dims,id_RUG)
232  call nf95_def_var(nid, "longitude", NF90_FLOAT, ndim, varid_longitude)
233  call nf95_def_var(nid, "latitude", NF90_FLOAT, ndim, varid_latitude)
234
235  !--- Attributes creation
236  ierr=NF90_PUT_ATT(nid,id_tim, "title","Jour dans l annee")
237  ierr=NF90_PUT_ATT(nid,id_FOCE,"title","Fraction ocean")
238  ierr=NF90_PUT_ATT(nid,id_FSIC,"title","Fraction glace de mer")
239  ierr=NF90_PUT_ATT(nid,id_FTER,"title","Fraction terre")
240  ierr=NF90_PUT_ATT(nid,id_FLIC,"title","Fraction land ice")
241  ierr=NF90_PUT_ATT(nid,id_SST ,"title","Temperature superficielle de la mer")
242  ierr=NF90_PUT_ATT(nid,id_BILS,"title","Reference flux de chaleur au sol")
243  ierr=NF90_PUT_ATT(nid,id_ALB, "title","Albedo a la surface")
244  ierr=NF90_PUT_ATT(nid,id_RUG, "title","Rugosite")
245
246  call nf95_put_att(nid, varid_longitude, "standard_name", "longitude")
247  call nf95_put_att(nid, varid_longitude, "units", "degrees_east")
248
249  call nf95_put_att(nid, varid_latitude, "standard_name", "latitude")
250  call nf95_put_att(nid, varid_latitude, "units", "degrees_north")
251
252  ierr=NF90_ENDDEF(nid)
253
254  !--- Variables saving
255  ierr=NF90_PUT_VAR(nid,id_tim,(/(REAL(k),k=1,ndays)/))
256  ierr=NF90_PUT_VAR(nid,id_FOCE,pctsrf_t(:,is_oce,:),(/1,1/),(/klon,ndays/))
257  ierr=NF90_PUT_VAR(nid,id_FSIC,pctsrf_t(:,is_sic,:),(/1,1/),(/klon,ndays/))
258  ierr=NF90_PUT_VAR(nid,id_FTER,pctsrf_t(:,is_ter,:),(/1,1/),(/klon,ndays/))
259  ierr=NF90_PUT_VAR(nid,id_FLIC,pctsrf_t(:,is_lic,:),(/1,1/),(/klon,ndays/))
260  ierr=NF90_PUT_VAR(nid,id_SST ,phy_sst(:,:),(/1,1/),(/klon,ndays/))
261  ierr=NF90_PUT_VAR(nid,id_BILS,phy_bil(:,:),(/1,1/),(/klon,ndays/))
262  ierr=NF90_PUT_VAR(nid,id_ALB ,phy_alb(:,:),(/1,1/),(/klon,ndays/))
263  ierr=NF90_PUT_VAR(nid,id_RUG ,phy_rug(:,:),(/1,1/),(/klon,ndays/))
264  call nf95_put_var(nid, varid_longitude, rlon)
265  call nf95_put_var(nid, varid_latitude, rlat)
266
267  ierr=NF90_CLOSE(nid)
268
269  IF (prt_level>5) WRITE(lunout,*) 'Ecriture du fichier limit : fin'
270
271  DEALLOCATE(pctsrf_t,phy_sst,phy_bil,phy_alb,phy_rug)
272
273
274!===============================================================================
275!
276  CONTAINS
277!
278!===============================================================================
279
280
281!-------------------------------------------------------------------------------
282!
283SUBROUTINE get_2Dfield(fnam, varname, mode, ibar, ndays, champo, flag, mask)
284!
285!-----------------------------------------------------------------------------
286! Comments:
287!   There are two assumptions concerning the NetCDF files, that are satisfied
288!   with files that are conforming NC convention:
289!     1) The last dimension of the variables used is the time record.
290!     2) Dimensional variables have the same names as corresponding dimensions.
291!-----------------------------------------------------------------------------
292  USE netcdf, ONLY: NF90_OPEN, NF90_INQ_VARID, NF90_INQUIRE_VARIABLE, &
293       NF90_CLOSE, NF90_INQ_DIMID, NF90_INQUIRE_DIMENSION, NF90_GET_VAR, &
294       NF90_GET_ATT
295  USE dimphy, ONLY : klon
296  USE phys_state_var_mod, ONLY : pctsrf
297  USE control_mod
298  USE pchsp_95_m, only: pchsp_95
299  USE pchfe_95_m, only: pchfe_95
300  USE arth_m, only: arth
301  USE indice_sol_mod
302
303  IMPLICIT NONE
304#include "dimensions.h"
305#include "paramet.h"
306#include "comgeom2.h"
307#include "iniprint.h"
308!-----------------------------------------------------------------------------
309! Arguments:
310  CHARACTER(LEN=*),  INTENT(IN)     :: fnam     ! NetCDF file name
311  CHARACTER(LEN=10), INTENT(IN)     :: varname  ! NetCDF variable name
312  CHARACTER(LEN=3),  INTENT(IN)     :: mode     ! RUG, SIC, SST or ALB
313  LOGICAL,           INTENT(IN)     :: ibar     ! interp on pressure levels
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,           DIMENSION(4) :: dids       ! NetCDF dimensions identifiers
325  REAL, ALLOCATABLE, DIMENSION(:) :: dlon_ini   ! initial longitudes vector
326  REAL, ALLOCATABLE, DIMENSION(:) :: dlat_ini   ! initial latitudes  vector
327  REAL, POINTER,     DIMENSION(:) :: dlon, dlat ! reordered lon/lat  vectors
328!--- fields
329  INTEGER :: imdep, jmdep, lmdep                ! dimensions of 'champ'
330  REAL, ALLOCATABLE, DIMENSION(:, :) :: champ   ! wanted field on initial grid
331  REAL, ALLOCATABLE, DIMENSION(:) :: yder, timeyear
332  REAL,              DIMENSION(iim, jjp1) :: champint   ! interpolated field
333  REAL, ALLOCATABLE, DIMENSION(:, :, :)    :: champtime
334  REAL, ALLOCATABLE, DIMENSION(:, :, :)    :: 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, DIMENSION(:, :) :: work     ! used for extrapolation
342  CHARACTER(LEN=25)                 :: title    ! for messages
343  LOGICAL                           :: extrp    ! flag for extrapolation
344  LOGICAL                           :: oldice   ! flag for old way ice computation
345  REAL                              :: chmin, chmax
346  INTEGER ierr
347  integer n_extrap ! number of extrapolated points
348  logical skip
349
350!------------------------------------------------------------------------------
351!---Variables depending on keyword 'mode' -------------------------------------
352  NULLIFY(champo)
353
354  SELECT CASE(mode)
355  CASE('RUG'); title='Rugosite'
356  CASE('SIC'); title='Sea-ice'
357  CASE('SST'); title='SST'
358  CASE('ALB'); title='Albedo'
359  END SELECT
360 
361
362  extrp=.FALSE.
363  oldice=.FALSE.
364  IF ( PRESENT(flag) ) THEN
365    IF ( flag .AND. mode=='SST' ) extrp=.TRUE.
366    IF ( flag .AND. mode=='SIC' ) oldice=.TRUE.
367  END IF
368
369!--- GETTING SOME DIMENSIONAL VARIABLES FROM FILE -----------------------------
370  IF (prt_level>5) WRITE(lunout,*) ' Now reading file : ',fnam
371  ierr=NF90_OPEN(fnam, NF90_NOWRITE, ncid);             CALL ncerr(ierr, fnam)
372  ierr=NF90_INQ_VARID(ncid, trim(varname), varid);            CALL ncerr(ierr, fnam)
373  ierr=NF90_INQUIRE_VARIABLE(ncid, varid, dimids=dids); CALL ncerr(ierr, fnam)
374
375!--- Read unit for sea-ice variable only
376  IF (mode=='SIC') THEN
377     ierr=NF90_GET_ATT(ncid, varid, 'units', unit_sic)
378     IF(ierr/=NF90_NOERR) THEN
379        IF (prt_level>5) WRITE(lunout,*) 'No unit was given in sea-ice file. Take percentage as default value'
380        unit_sic='X'
381     ELSE
382        IF (prt_level>5) WRITE(lunout,*) ' Sea-ice cover has unit=',unit_sic
383     END IF
384  END IF
385
386!--- Longitude
387  ierr=NF90_INQUIRE_DIMENSION(ncid, dids(1), name=dnam, len=imdep)
388  CALL ncerr(ierr, fnam); ALLOCATE(dlon_ini(imdep), dlon(imdep))
389  ierr=NF90_INQ_VARID(ncid, dnam, varid);                CALL ncerr(ierr, fnam)
390  ierr=NF90_GET_VAR(ncid, varid, dlon_ini);              CALL ncerr(ierr, fnam)
391  IF (prt_level>5) WRITE(lunout, *) 'variable ', dnam, 'dimension ', imdep
392
393!--- Latitude
394  ierr=NF90_INQUIRE_DIMENSION(ncid, dids(2), name=dnam, len=jmdep)
395  CALL ncerr(ierr, fnam); ALLOCATE(dlat_ini(jmdep), dlat(jmdep))
396  ierr=NF90_INQ_VARID(ncid, dnam, varid);                CALL ncerr(ierr, fnam)
397  ierr=NF90_GET_VAR(ncid, varid, dlat_ini);              CALL ncerr(ierr, fnam)
398  IF (prt_level>5) WRITE(lunout, *) 'variable ', dnam, 'dimension ', jmdep
399
400!--- Time (variable is not needed - it is rebuilt - but calendar is)
401  ierr=NF90_INQUIRE_DIMENSION(ncid, dids(3), name=dnam, len=lmdep)
402  CALL ncerr(ierr, fnam); ALLOCATE(timeyear(lmdep))
403  ierr=NF90_INQ_VARID(ncid, dnam, varid);                CALL ncerr(ierr, fnam)
404  cal_in=' '
405  ierr=NF90_GET_ATT(ncid, varid, 'calendar', cal_in)
406  IF(ierr/=NF90_NOERR) THEN
407    SELECT CASE(mode)
408      CASE('RUG', 'ALB'); cal_in='360d'
409      CASE('SIC', 'SST'); cal_in='gregorian'
410    END SELECT
411    IF (prt_level>5) WRITE(lunout, *)'ATTENTION: variable "time" sans attribut "calendrier" ' &
412         // 'dans '//TRIM(fnam)//'. On choisit la valeur par defaut.'
413  END IF
414  IF (prt_level>5) WRITE(lunout, *) 'variable ', dnam, 'dimension ', lmdep, 'calendrier ', &
415       cal_in
416
417 
418!--- CONSTRUCTING THE INPUT TIME VECTOR FOR INTERPOLATION --------------------
419  !--- Determining input file number of days, depending on calendar
420  ndays_in=year_len(anneeref, cal_in)
421
422!--- Time vector reconstruction (time vector from file is not trusted)
423!--- If input records are not monthly, time sampling has to be constant !
424  timeyear=mid_months(anneeref, cal_in, lmdep)
425  IF (lmdep /= 12) WRITE(lunout,*) 'Note : les fichiers de ', TRIM(mode), &
426       ' ne comportent pas 12, mais ', lmdep, ' enregistrements.'
427
428!--- GETTING THE FIELD AND INTERPOLATING IT ----------------------------------
429  ALLOCATE(champ(imdep, jmdep), champtime(iim, jjp1, lmdep))
430  IF(extrp) ALLOCATE(work(imdep, jmdep))
431
432  IF (prt_level>5) WRITE(lunout, *)
433  IF (prt_level>5) WRITE(lunout,*)'LECTURE ET INTERPOLATION HORIZ. DE ', lmdep, ' CHAMPS.'
434  ierr=NF90_INQ_VARID(ncid, varname, varid);             CALL ncerr(ierr, fnam)
435  DO l=1, lmdep
436    ierr=NF90_GET_VAR(ncid, varid, champ, (/1, 1, l/), (/imdep, jmdep, 1/))
437    CALL ncerr(ierr, fnam)
438    CALL conf_dat2d(title, imdep, jmdep, dlon_ini, dlat_ini, dlon, dlat, &
439         champ, ibar)
440
441    IF (extrp) CALL extrapol(champ, imdep, jmdep, 999999., .TRUE., .TRUE., 2, &
442         work)
443
444    IF(ibar .AND. .NOT.oldice) THEN
445      IF(l==1 .AND. prt_level>5) THEN
446        WRITE(lunout, *) '-------------------------------------------------------------------------'
447        WRITE(lunout, *) 'Utilisation de l''interpolation barycentrique pour ',TRIM(title),' $$$'
448        WRITE(lunout, *) '-------------------------------------------------------------------------'
449      END IF
450      IF(mode=='RUG') champ=LOG(champ)
451      CALL inter_barxy(dlon, dlat(:jmdep-1), champ, rlonu(:iim),     &
452                         rlatv, champint)
453      IF(mode=='RUG') THEN
454        champint=EXP(champint)
455        WHERE(NINT(mask)/=1) champint=0.001
456      END IF
457    ELSE
458      SELECT CASE(mode)
459        CASE('RUG');       CALL rugosite(imdep, jmdep, dlon, dlat, champ,    &
460                                    iim, jjp1, rlonv, rlatu, champint, mask)
461        CASE('SIC');       CALL sea_ice (imdep, jmdep, dlon, dlat, champ,    &
462                                    iim, jjp1, rlonv, rlatu, champint)
463        CASE('SST', 'ALB'); CALL grille_m(imdep, jmdep, dlon, dlat, champ,    &
464                                    iim, jjp1, rlonv, rlatu, champint)
465      END SELECT
466    END IF
467    champtime(:, :, l)=champint
468  END DO
469  ierr=NF90_CLOSE(ncid);                                 CALL ncerr(ierr, fnam)
470
471  DEALLOCATE(dlon_ini, dlat_ini, dlon, dlat, champ)
472  IF(extrp) DEALLOCATE(work)
473
474!--- TIME INTERPOLATION ------------------------------------------------------
475  IF (prt_level>5) THEN
476     WRITE(lunout, *)
477     WRITE(lunout, *)'INTERPOLATION TEMPORELLE.'
478     WRITE(lunout, *)' Vecteur temps en entree: ', timeyear
479     WRITE(lunout, *)' Vecteur temps en sortie de 0 a ', ndays
480  END IF
481
482  ALLOCATE(yder(lmdep), champan(iip1, jjp1, ndays))
483  skip = .false.
484  n_extrap = 0
485  DO j=1, jjp1
486    DO i=1, iim
487      yder = pchsp_95(timeyear, champtime(i, j, :), ibeg=2, iend=2, &
488           vc_beg=0., vc_end=0.)
489      CALL pchfe_95(timeyear, champtime(i, j, :), yder, skip, &
490           arth(0., real(ndays_in) / ndays, ndays), champan(i, j, :), ierr)
491      if (ierr < 0) stop 1
492      n_extrap = n_extrap + ierr
493    END DO
494  END DO
495  if (n_extrap /= 0) then
496     WRITE(lunout,*) "get_2Dfield pchfe_95: n_extrap = ", n_extrap
497  end if
498  champan(iip1, :, :)=champan(1, :, :)
499  DEALLOCATE(yder, champtime, timeyear)
500
501!--- Checking the result
502  DO j=1, jjp1
503    CALL minmax(iip1, champan(1, j, 10), chmin, chmax)
504    IF (prt_level>5) WRITE(lunout, *)' ',TRIM(title),' au temps 10 ', chmin, chmax, j
505  END DO
506
507!--- SPECIAL FILTER FOR SST: SST>271.38 --------------------------------------
508  IF(mode=='SST') THEN
509    IF (prt_level>5) WRITE(lunout, *) 'Filtrage de la SST: SST >= 271.38'
510    WHERE(champan<271.38) champan=271.38
511  END IF
512
513!--- SPECIAL FILTER FOR SIC: 0.0<SIC<1.0 -------------------------------------
514  IF(mode=='SIC') THEN
515    IF (prt_level>5) WRITE(lunout, *) 'Filtrage de la SIC: 0.0 < Sea-ice < 1.0'
516
517    IF (unit_sic=='1') THEN
518       ! Nothing to be done for sea-ice field is already in fraction of 1
519       ! This is the case for sea-ice in file cpl_atm_sic.nc
520       IF (prt_level>5) WRITE(lunout,*) 'Sea-ice field already in fraction of 1'
521    ELSE
522       ! Convert sea ice from percentage to fraction of 1
523       IF (prt_level>5) WRITE(lunout,*) 'Transformt sea-ice field from percentage to fraction of 1.'
524       champan(:, :, :)=champan(:, :, :)/100.
525    END IF
526
527    champan(iip1, :, :)=champan(1, :, :)
528    WHERE(champan>1.0) champan=1.0
529    WHERE(champan<0.0) champan=0.0
530 END IF
531
532!--- DYNAMICAL TO PHYSICAL GRID ----------------------------------------------
533  ALLOCATE(champo(klon, ndays))
534  DO k=1, ndays
535    CALL gr_dyn_fi(1, iip1, jjp1, klon, champan(1, 1, k), champo(1, k))
536  END DO
537  DEALLOCATE(champan)
538
539END SUBROUTINE get_2Dfield
540!
541!-------------------------------------------------------------------------------
542
543
544
545!-------------------------------------------------------------------------------
546!
547FUNCTION year_len(y,cal_in)
548!
549!-------------------------------------------------------------------------------
550  USE ioipsl, ONLY : ioget_calendar,ioconf_calendar,lock_calendar,ioget_year_len
551  IMPLICIT NONE
552!-------------------------------------------------------------------------------
553! Arguments:
554  INTEGER                       :: year_len
555  INTEGER,           INTENT(IN) :: y
556  CHARACTER(LEN=*),  INTENT(IN) :: cal_in
557!-------------------------------------------------------------------------------
558! Local variables:
559  CHARACTER(LEN=20)             :: cal_out              ! calendar (for outputs)
560!-------------------------------------------------------------------------------
561!--- Getting the input calendar to reset at the end of the function
562  CALL ioget_calendar(cal_out)
563
564!--- Unlocking calendar and setting it to wanted one
565  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_in))
566
567!--- Getting the number of days in this year
568  year_len=ioget_year_len(y)
569
570!--- Back to original calendar
571  CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_out))
572
573END FUNCTION year_len
574!
575!-------------------------------------------------------------------------------
576
577
578!-------------------------------------------------------------------------------
579!
580FUNCTION mid_months(y,cal_in,nm)
581!
582!-------------------------------------------------------------------------------
583  USE ioipsl, ONLY : ioget_calendar,ioconf_calendar,lock_calendar,ioget_mon_len
584  IMPLICIT NONE
585!-------------------------------------------------------------------------------
586! Arguments:
587  INTEGER,                INTENT(IN) :: y               ! year
588  CHARACTER(LEN=*),       INTENT(IN) :: cal_in          ! calendar
589  INTEGER,                INTENT(IN) :: nm              ! months/year number
590  REAL,    DIMENSION(nm)             :: mid_months      ! mid-month times
591!-------------------------------------------------------------------------------
592! Local variables:
593  CHARACTER(LEN=99)                  :: mess            ! error message
594  CHARACTER(LEN=20)                  :: cal_out         ! calendar (for outputs)
595  INTEGER, DIMENSION(nm)             :: mnth            ! months lengths (days)
596  INTEGER                            :: m               ! months counter
597  INTEGER                            :: nd              ! number of days
598!-------------------------------------------------------------------------------
599  nd=year_len(y,cal_in)
600
601  IF(nm==12) THEN
602
603  !--- Getting the input calendar to reset at the end of the function
604    CALL ioget_calendar(cal_out)
605
606  !--- Unlocking calendar and setting it to wanted one
607    CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_in))
608
609  !--- Getting the length of each month
610    DO m=1,nm; mnth(m)=ioget_mon_len(y,m); END DO
611
612  !--- Back to original calendar
613    CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_out))
614
615  ELSE IF(MODULO(nd,nm)/=0) THEN
616    WRITE(mess,'(a,i3,a,i3,a)')'Unconsistent calendar: ',nd,' days/year, but ',&
617      nm,' months/year. Months number should divide days number.'
618    CALL abort_gcm('mid_months',TRIM(mess),1)
619
620  ELSE
621    mnth=(/(m,m=1,nm,nd/nm)/)
622  END IF
623
624!--- Mid-months times
625  mid_months(1)=0.5*REAL(mnth(1))
626  DO k=2,nm
627    mid_months(k)=mid_months(k-1)+0.5*REAL(mnth(k-1)+mnth(k))
628  END DO
629
630END FUNCTION mid_months
631!
632!-------------------------------------------------------------------------------
633
634
635!-------------------------------------------------------------------------------
636!
637SUBROUTINE ncerr(ncres,fnam)
638!
639!-------------------------------------------------------------------------------
640! Purpose: NetCDF errors handling.
641!-------------------------------------------------------------------------------
642  USE netcdf, ONLY : NF90_NOERR, NF90_STRERROR
643  IMPLICIT NONE
644!-------------------------------------------------------------------------------
645! Arguments:
646  INTEGER,          INTENT(IN) :: ncres
647  CHARACTER(LEN=*), INTENT(IN) :: fnam
648!-------------------------------------------------------------------------------
649#include "iniprint.h"
650  IF(ncres/=NF90_NOERR) THEN
651    WRITE(lunout,*)'Problem with file '//TRIM(fnam)//' in routine limit_netcdf.'
652    CALL abort_gcm('limit_netcdf',NF90_STRERROR(ncres),1)
653  END IF
654
655END SUBROUTINE ncerr
656!
657!-------------------------------------------------------------------------------
658
659#endif
660! of #ifdef CPP_EARTH
661
662END SUBROUTINE limit_netcdf
Note: See TracBrowser for help on using the repository browser.