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

Last change on this file since 2146 was 2128, checked in by lguez, 10 years ago

Added longitude and latitude to file "limit.nc", so that we can more easily
plot fields from this file.

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