source: LMDZ5/trunk/libf/dynlonlat_phylonlat/phylmd/etat0dyn_netcdf.F90 @ 2354

Last change on this file since 2354 was 2336, checked in by dcugnet, 9 years ago
  • Add parallel capability for ce0l.
  • Small bug in grid_noro fixed (smoothed topography was used instead of unsmoothed one for geopotential computation at north pole).
  • Removed average of mass at poles in etat0dyn_netcdf after start_init_dyn => different results in the zoomed grid case.
  • ok_etat0=n and ok_limit=y combination now works fine (if no initial state is needed, but only limit.nc file). This required:
    • to move grid_noro0 and start_init_noro0 subroutines from etat0dyn_netcdf.F90 to limit_netcdf.F90
    • to create init_ssrf_m.F90 file, so that sub-surfaces can be initialized from limit_netcdf.F90 without any etat0*_netcdf routines call).
  • Simplified somehow the corresponding code, in particular: 1) removed obsolete flags "oldice". 2) removed flag "ibar": barycentric interpolation is used everywhere (except in start_init_subsurf, still calling grille_m - to be changed soon). 3) removed useless CPP_PHY precompilation directives, considering the possibility to run ce0l without physics is useless (ce0l is dedicated to Earth physics).
File size: 22.2 KB
Line 
1MODULE etat0dyn
2!
3!*******************************************************************************
4! Purpose: Create dynamical initial state using atmospheric fields from a
5!          database of atmospheric to initialize the model.
6!-------------------------------------------------------------------------------
7! Comments:
8!
9!    *  This module is designed to work for Earth (and with ioipsl)
10!
11!    *  etat0dyn_netcdf routine can access to NetCDF data through the following
12!  routine (to be called after restget):
13!    CALL startget_dyn3d(varname, lon_in,  lat_in, pls, workvar,&
14!                          champ, lon_in2, lat_in2)
15!
16!    *  Variables should have the following names in the NetCDF files:
17!            'U'      : East ward wind              (in "ECDYN.nc")
18!            'V'      : Northward wind              (in "ECDYN.nc")
19!            'TEMP'   : Temperature                 (in "ECDYN.nc")
20!            'R'      : Relative humidity           (in "ECDYN.nc")
21!            'RELIEF' : High resolution orography   (in "Relief.nc")
22!
23!    * The land mask and corresponding weights can be:
24!      1) already known (in particular if etat0dyn has been called before) ;
25!         in this case, ANY(masque(:,:)/=-99999.) = .TRUE.
26!      2) computed using the ocean mask from the ocean model (to ensure ocean
27!         fractions are the same for atmosphere and ocean) for coupled runs.
28!         File name: "o2a.nc"  ;  variable name: "OceMask"
29!      3) computed from topography file "Relief.nc" for forced runs.
30!
31!   *   There is a big mess with the longitude size. Should it be iml or iml+1 ?
32!  I have chosen to use the iml+1 as an argument to this routine and we declare
33!  internaly smaller fields when needed. This needs to be cleared once and for
34!  all in LMDZ. A convention is required.
35!-------------------------------------------------------------------------------
36  USE ioipsl,         ONLY: flininfo, flinopen, flinget, flinclo, histclo
37  USE assert_eq_m,    ONLY: assert_eq
38  USE indice_sol_mod, ONLY: epsfra
39  IMPLICIT NONE
40
41  PRIVATE
42  PUBLIC :: etat0dyn_netcdf
43
44  include "iniprint.h"
45  include "dimensions.h"
46  include "paramet.h"
47  include "comgeom2.h"
48  include "comvert.h"
49  include "comconst.h"
50  include "temps.h"
51  include "comdissnew.h"
52  include "serre.h"
53  REAL, SAVE :: deg2rad
54  INTEGER,            SAVE      :: iml_dyn, jml_dyn, llm_dyn, ttm_dyn, fid_dyn
55  REAL, ALLOCATABLE,  SAVE      :: lon_dyn(:,:), lat_dyn(:,:), levdyn_ini(:)
56  CHARACTER(LEN=120), PARAMETER :: dynfname='ECDYN.nc'
57
58CONTAINS
59
60!-------------------------------------------------------------------------------
61!
62SUBROUTINE etat0dyn_netcdf(masque, phis)
63!
64!-------------------------------------------------------------------------------
65! Purpose: Create dynamical initial states.
66!-------------------------------------------------------------------------------
67! Notes:  1) This routine is designed to work for Earth
68!         2) If masque(:,:)/=-99999., masque and phis are already known.
69!         Otherwise: compute it.
70!-------------------------------------------------------------------------------
71  USE control_mod
72  USE regr_lat_time_coefoz_m, ONLY: regr_lat_time_coefoz
73  USE regr_pr_o3_m,   ONLY: regr_pr_o3
74  USE press_coefoz_m, ONLY: press_coefoz
75  USE exner_hyb_m,    ONLY: exner_hyb
76  USE exner_milieu_m, ONLY: exner_milieu
77  USE infotrac,       ONLY: nqtot, tname
78  USE filtreg_mod
79  IMPLICIT NONE
80!-------------------------------------------------------------------------------
81! Arguments:
82  REAL,    INTENT(INOUT) :: masque(iip1,jjp1)   !--- Land-ocean mask
83  REAL,    INTENT(INOUT) :: phis  (iip1,jjp1)   !--- Ground geopotential
84!-------------------------------------------------------------------------------
85! Local variables:
86  CHARACTER(LEN=256) :: modname, fmt
87  INTEGER            :: i, j, l, ji, itau, iday
88  REAL               :: xpn, xps, time, phystep
89  REAL, DIMENSION(iip1,jjp1)       :: psol, masque_tmp
90  REAL, DIMENSION(iip1,jjp1,llm+1) :: p3d
91  REAL, DIMENSION(iip1,jjp1,llm)   :: uvent, t3d, tpot, qsat, qd
92  REAL, DIMENSION(iip1,jjp1,llm)   :: pk, pls, y, masse
93  REAL, DIMENSION(iip1,jjm ,llm)   :: vvent
94  REAL, DIMENSION(ip1jm    ,llm)   :: pbarv
95  REAL, DIMENSION(ip1jmp1  ,llm)   :: pbaru, phi, w
96  REAL, DIMENSION(ip1jmp1)         :: pks
97  REAL, DIMENSION(iim)             :: xppn, xpps
98  REAL, ALLOCATABLE                :: q3d(:,:,:,:)
99!-------------------------------------------------------------------------------
100  modname='etat0dyn_netcdf'
101
102  deg2rad = pi/180.0
103
104! Compute ground geopotential and possibly the mask.
105!*******************************************************************************
106  masque_tmp(:,:)=masque(:,:)
107  WRITE(fmt,"(i4,'i1)')")iip1 ; fmt='('//ADJUSTL(fmt)
108  IF(ALL(masque==-99999.)) THEN                         !--- KEEP NEW MASK
109    masque=masque_tmp
110    IF(prt_level>=1) THEN
111      WRITE(lunout,*)'BUILT MASK :'
112      WRITE(lunout,fmt) NINT(masque)
113    END IF
114    WHERE(   masque(:,:)<EPSFRA) masque(:,:)=0.
115    WHERE(1.-masque(:,:)<EPSFRA) masque(:,:)=1.
116  END IF
117
118! Compute psol AND tsol, knowing phis.
119!*******************************************************************************
120  CALL start_init_dyn(rlonv, rlatu, rlonu, rlatv, phis, psol)
121
122! Mid-levels pressure computation
123!*******************************************************************************
124  CALL pression(ip1jmp1, ap, bp, psol, p3d)             !--- Update p3d
125  IF(pressure_exner) THEN                               !--- Update pk, pks
126    CALL exner_hyb   (ip1jmp1,psol,p3d,pks,pk)
127  ELSE
128    CALL exner_milieu(ip1jmp1,psol,p3d,pks,pk)
129  END IF
130  pls(:,:,:)=preff*(pk(:,:,:)/cpp)**(1./kappa)          !--- Update pls
131
132! Update uvent, vvent, t3d and tpot
133!*******************************************************************************
134  uvent(:,:,:) = 0.0 ; vvent(:,:,:) = 0.0 ; t3d (:,:,:) = 0.0
135  CALL startget_dyn3d('u'   ,rlonu,rlatu,pls,y ,uvent,rlonv,rlatv)
136  CALL startget_dyn3d('v'   ,rlonv,rlatv,pls(:,:jjm,:),y(:,:jjm,:),vvent,      &
137 &                           rlonu,rlatu(:jjm))
138  CALL startget_dyn3d('t'   ,rlonv,rlatu,pls,y ,t3d ,rlonu,rlatv)
139  tpot(:,:,:)=t3d(:,:,:)
140  CALL startget_dyn3d('tpot',rlonv,rlatu,pls,pk,tpot,rlonu,rlatv)
141
142  WRITE(lunout,*) 'T3D min,max:',MINVAL(t3d(:,:,:)),MAXVAL(t3d(:,:,:))
143  WRITE(lunout,*) 'PLS min,max:',MINVAL(pls(:,:,:)),MAXVAL(pls(:,:,:))
144
145! Humidity at saturation computation
146!*******************************************************************************
147  WRITE(lunout,*) 'avant q_sat'
148  CALL q_sat(llm*jjp1*iip1, t3d, pls, qsat)
149  WRITE(lunout,*) 'apres q_sat'
150  WRITE(lunout,*) 'QSAT min,max:',MINVAL(qsat(:,:,:)),MAXVAL(qsat(:,:,:))
151!  WRITE(lunout,*) 'QSAT :',qsat(10,20,:)
152  qd (:,:,:) = 0.0
153  CALL startget_dyn3d('q',rlonv,rlatu,pls,qsat,qd,rlonu,rlatv)
154  ALLOCATE(q3d(iip1,jjp1,llm,nqtot)); q3d(:,:,:,:)=0.0 ; q3d(:,:,:,1)=qd(:,:,:)
155  CALL flinclo(fid_dyn)
156
157! Parameterization of ozone chemistry:
158!*******************************************************************************
159! Look for ozone tracer:
160  DO i=1,nqtot; IF(ANY(["O3","o3"]==tname(i))) EXIT; END DO
161  IF(i/=nqtot+1) THEN
162    CALL regr_lat_time_coefoz
163    CALL press_coefoz
164    CALL regr_pr_o3(p3d, q3d(:,:,:,i))
165    q3d(:,:,:,i)=q3d(:,:,:,i)*48./ 29.                  !--- Mole->mass fraction         
166  END IF
167  q3d(iip1,:,:,:)=q3d(1,:,:,:)
168
169! Writing
170!*******************************************************************************
171  CALL inidissip(lstardis, nitergdiv, nitergrot, niterh, tetagdiv, tetagrot,   &
172       tetatemp, vert_prof_dissip)
173  WRITE(lunout,*)'sortie inidissip'
174  itau=0
175  itau_dyn=0
176  itau_phy=0
177  iday=dayref+itau/day_step
178  time=FLOAT(itau-(iday-dayref)*day_step)/day_step
179  IF(time>1.) THEN
180   time=time-1
181   iday=iday+1
182  END IF
183  day_ref=dayref
184  annee_ref=anneeref
185  CALL geopot( ip1jmp1, tpot, pk, pks, phis, phi )
186  WRITE(lunout,*)'sortie geopot'
187  CALL caldyn0( itau, uvent, vvent, tpot, psol, masse, pk, phis,               &
188                phi,  w, pbaru, pbarv, time+iday-dayref)
189  WRITE(lunout,*)'sortie caldyn0'     
190#ifdef CPP_PARA
191  CALL dynredem0_loc( "start.nc", dayref, phis)
192#else
193  CALL dynredem0( "start.nc", dayref, phis)
194#endif
195  WRITE(lunout,*)'sortie dynredem0'
196#ifdef CPP_PARA
197  CALL dynredem1_loc( "start.nc", 0.0, vvent, uvent, tpot, q3d, masse, psol)
198#else
199  CALL dynredem1( "start.nc", 0.0, vvent, uvent, tpot, q3d, masse, psol)
200#endif
201  WRITE(lunout,*)'sortie dynredem1'
202  CALL histclo()
203
204END SUBROUTINE etat0dyn_netcdf
205!
206!-------------------------------------------------------------------------------
207
208
209!-------------------------------------------------------------------------------
210!
211SUBROUTINE startget_dyn3d(var, lon_in,  lat_in,  pls,  workvar,&
212                        champ, lon_in2, lat_in2)
213!-------------------------------------------------------------------------------
214  IMPLICIT NONE
215!===============================================================================
216! Purpose: Compute some quantities (u,v,t,q,tpot) using variables U,V,TEMP and R
217!     (3D fields) of file dynfname.
218!-------------------------------------------------------------------------------
219! Note: An input auxilliary field "workvar" has to be specified in two cases:
220!     * for "q":    the saturated humidity.
221!     * for "tpot": the Exner function.
222!===============================================================================
223! Arguments:
224  CHARACTER(LEN=*), INTENT(IN)    :: var
225  REAL,             INTENT(IN)    :: lon_in(:)        ! dim (iml)
226  REAL,             INTENT(IN)    :: lat_in(:)        ! dim (jml)
227  REAL,             INTENT(IN)    :: pls    (:, :, :) ! dim (iml, jml, lml)
228  REAL,             INTENT(IN)    :: workvar(:, :, :) ! dim (iml, jml, lml)
229  REAL,             INTENT(INOUT) :: champ  (:, :, :) ! dim (iml, jml, lml)
230  REAL,             INTENT(IN)    :: lon_in2(:)       ! dim (iml)
231  REAL,             INTENT(IN)    :: lat_in2(:)       ! dim (jml2)
232!-------------------------------------------------------------------------------
233! Local variables:
234  CHARACTER(LEN=10)  :: vname
235  CHARACTER(LEN=256) :: msg, modname="startget_dyn3d"
236  INTEGER            :: iml, jml, jml2, lml, il
237  REAL               :: xppn, xpps
238!-------------------------------------------------------------------------------
239  iml=assert_eq([SIZE(lon_in),SIZE(pls,1),SIZE(workvar,1),SIZE(champ,1),       &
240     &                                    SIZE(lon_in2)], TRIM(modname)//" iml")
241  jml=assert_eq( SIZE(lat_in),SIZE(pls,2),SIZE(workvar,2),SIZE(champ,2),       &
242     &                                                    TRIM(modname)//" jml")
243  lml=assert_eq(              SIZE(pls,3),SIZE(workvar,3),SIZE(champ,3),       &
244     &                                                    TRIM(modname)//" lml")
245  jml2=SIZE(lat_in2)
246
247!--- CHECK IF THE FIELD IS KNOWN
248   SELECT CASE(var)
249    CASE('u');    vname='U'
250    CASE('v');    vname='V'
251    CASE('t');    vname='TEMP'
252    CASE('q');    vname='R';    msg='humidity as the saturated humidity'
253    CASE('tpot'); msg='potential temperature as the Exner function'
254    CASE DEFAULT; msg='No rule to extract variable '//TRIM(var)
255      CALL abort_gcm(modname,TRIM(msg)//' from any data set',1)
256  END SELECT
257
258!--- CHECK IF SOMETHING IS MISSING
259  IF((var=='tpot'.OR.var=='q').AND.MINVAL(workvar)==MAXVAL(workvar)) THEN
260    msg='Could not compute '//TRIM(msg)//' is missing or constant.'
261    CALL abort_gcm(modname,TRIM(msg),1)
262  END IF
263
264!--- INTERPOLATE 3D FIELD IF NEEDED
265  IF(var/='tpot') CALL start_inter_3d(TRIM(vname),lon_in,lat_in,lon_in2,      &
266                                                  lat_in2,pls,champ)
267
268!--- COMPUTE THE REQUIRED FILED
269  SELECT CASE(var)
270    CASE('u'); DO il=1,lml; champ(:,:,il)=champ(:,:,il)*cu(:,1:jml); END DO
271      champ(iml,:,:)=champ(1,:,:)                   !--- Eastward wind
272
273    CASE('v'); DO il=1,lml; champ(:,:,il)=champ(:,:,il)*cv(:,1:jml); END DO
274      champ(iml,:,:)=champ(1,:,:)                   !--- Northward wind
275
276    CASE('tpot','q')
277      IF(var=='tpot') THEN; champ=champ*cpp/workvar !--- Potential temperature
278      ELSE;                 champ=champ*.01*workvar !--- Relative humidity
279        WHERE(champ<0.) champ=1.0E-10
280      END IF
281      DO il=1,lml
282        xppn = SUM(aire(:,1  )*champ(:,1  ,il))/apoln
283        xpps = SUM(aire(:,jml)*champ(:,jml,il))/apols
284        champ(:,1  ,il) = xppn
285        champ(:,jml,il) = xpps
286      END DO
287  END SELECT
288
289END SUBROUTINE startget_dyn3d
290!
291!-------------------------------------------------------------------------------
292
293
294!-------------------------------------------------------------------------------
295!
296SUBROUTINE start_init_dyn(lon_in,lat_in,lon_in2,lat_in2,zs,psol)
297!
298!-------------------------------------------------------------------------------
299  IMPLICIT NONE
300!===============================================================================
301! Purpose:   Compute psol, knowing phis.
302!===============================================================================
303! Arguments:
304  REAL,    INTENT(IN)  :: lon_in (:),  lat_in (:)    ! dim (iml) (jml)
305  REAL,    INTENT(IN)  :: lon_in2(:),  lat_in2(:)    ! dim (iml) (jml2)
306  REAL,    INTENT(IN)  :: zs  (:,:)                  ! dim (iml,jml)
307  REAL,    INTENT(OUT) :: psol(:,:)                  ! dim (iml,jml)
308!-------------------------------------------------------------------------------
309! Local variables:
310  CHARACTER(LEN=256) :: modname='start_init_dyn'
311  REAL               :: date, dt
312  INTEGER            :: iml, jml, jml2, itau(1)
313  REAL, ALLOCATABLE  :: lon_rad(:), lon_ini(:), var_ana(:,:)
314  REAL, ALLOCATABLE  :: lat_rad(:), lat_ini(:)
315  REAL, ALLOCATABLE  :: z(:,:), ps(:,:), ts(:,:)
316!-------------------------------------------------------------------------------
317  iml=assert_eq(SIZE(lon_in),SIZE(zs,1),SIZE(psol,1),SIZE(lon_in2),            &
318      &                                              TRIM(modname)//" iml")
319  jml=assert_eq(SIZE(lat_in),SIZE(zs,2),SIZE(psol,2),TRIM(modname)//" jml")
320  jml2=SIZE(lat_in2)
321
322  WRITE(lunout,*) 'Opening the surface analysis'
323  CALL flininfo(dynfname, iml_dyn, jml_dyn, llm_dyn, ttm_dyn, fid_dyn)
324  WRITE(lunout,*) 'Values read: ', iml_dyn, jml_dyn, llm_dyn, ttm_dyn
325
326  ALLOCATE(lon_dyn(iml_dyn,jml_dyn), lat_dyn(iml_dyn,jml_dyn))
327  ALLOCATE(levdyn_ini(llm_dyn))
328  CALL flinopen(dynfname, .FALSE., iml_dyn, jml_dyn, llm_dyn,                  &
329                lon_dyn,lat_dyn,levdyn_ini,ttm_dyn,itau,date,dt,fid_dyn)
330
331!--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS
332  ALLOCATE(lon_ini(iml_dyn),lat_ini(jml_dyn))
333  lon_ini(:)=lon_dyn(:,1); IF(MAXVAL(lon_dyn)>pi) lon_ini=lon_ini*deg2rad
334  lat_ini(:)=lat_dyn(1,:); IF(MAXVAL(lat_dyn)>pi) lat_ini=lat_ini*deg2rad
335
336  ALLOCATE(var_ana(iml_dyn,jml_dyn),lon_rad(iml_dyn),lat_rad(jml_dyn))
337  CALL get_var_dyn('Z',z)                        !--- SURFACE GEOPOTENTIAL
338  CALL get_var_dyn('SP',ps)                      !--- SURFACE PRESSURE
339  CALL get_var_dyn('ST',ts)                      !--- SURFACE TEMPERATURE
340!  CALL flinclo(fid_dyn)
341  DEALLOCATE(var_ana,lon_rad,lat_rad,lon_ini,lat_ini)
342
343!--- PSOL IS COMPUTED IN PASCALS
344  psol(:iml-1,:) = ps(:iml-1,:)*(1.0+(z(:iml-1,:)-zs(:iml-1,:))/287.0          &
345     &            /ts(:iml-1,:))
346  psol(iml,:)=psol(1,:)
347  DEALLOCATE(z,ps,ts)
348  psol(:,1  )=SUM(aire(1:iml-1,1  )*psol(1:iml-1,1  ))/apoln  !--- NORTH POLE
349  psol(:,jml)=SUM(aire(1:iml-1,jml)*psol(1:iml-1,jml))/apols  !--- SOUTH POLE
350
351CONTAINS
352
353!-------------------------------------------------------------------------------
354!
355SUBROUTINE get_var_dyn(title,field)
356!
357!-------------------------------------------------------------------------------
358  USE conf_dat_m, ONLY: conf_dat2d
359  IMPLICIT NONE
360!-------------------------------------------------------------------------------
361! Arguments:
362  CHARACTER(LEN=*),  INTENT(IN)    :: title
363  REAL, ALLOCATABLE, INTENT(INOUT) :: field(:,:)
364!-------------------------------------------------------------------------------
365! Local variables:
366  CHARACTER(LEN=256) :: msg
367  INTEGER :: tllm
368!-------------------------------------------------------------------------------
369  SELECT CASE(title)
370    CASE('Z');     tllm=0;       msg='geopotential'
371    CASE('SP');    tllm=0;       msg='surface pressure'
372    CASE('ST');    tllm=llm_dyn; msg='temperature'
373  END SELECT
374  IF(.NOT.ALLOCATED(field)) THEN
375    ALLOCATE(field(iml,jml))
376    CALL flinget(fid_dyn, title, iml_dyn,jml_dyn, tllm, ttm_dyn, 1, 1, var_ana)
377    CALL conf_dat2d(title, lon_ini, lat_ini, lon_rad, lat_rad, var_ana, .TRUE.)
378    CALL interp_startvar(title, .TRUE., lon_rad,lat_rad, var_ana,              &
379                                        lon_in, lat_in, lon_in2, lat_in2, field)
380  ELSE IF(SIZE(field)/=SIZE(z)) THEN
381    msg='The '//TRIM(msg)//' field we have does not have the right size'
382    CALL abort_gcm(TRIM(modname),msg,1)
383  END IF
384
385END SUBROUTINE get_var_dyn
386!
387!-------------------------------------------------------------------------------
388
389END SUBROUTINE start_init_dyn
390!
391!-------------------------------------------------------------------------------
392
393
394!-------------------------------------------------------------------------------
395!
396SUBROUTINE start_inter_3d(var,lon_in,lat_in,lon_in2,lat_in2,pls_in,var3d)
397!
398!-------------------------------------------------------------------------------
399  USE conf_dat_m, ONLY: conf_dat3d
400  USE pchsp_95_m, ONLY: pchsp_95
401  USE pchfe_95_m, ONLY: pchfe_95
402  IMPLICIT NONE
403!-------------------------------------------------------------------------------
404! Arguments:
405  CHARACTER(LEN=*), INTENT(IN) :: var
406  REAL,    INTENT(IN)  :: lon_in(:),  lat_in(:)   ! dim (iml) (jml)
407  REAL,    INTENT(IN)  :: lon_in2(:), lat_in2(:)  ! dim (iml) (jml2)
408  REAL,    INTENT(IN)  :: pls_in(:,:,:)           ! dim (iml,jml,lml)
409  REAL,    INTENT(OUT) :: var3d (:,:,:)           ! dim (iml,jml,lml)
410!-------------------------------------------------------------------------------
411! Local variables:
412  CHARACTER(LEN=256) :: modname='start_inter_3d'
413  LOGICAL :: skip
414  REAL    :: chmin, chmax
415  INTEGER :: iml, jml, lml, jml2, ii, ij, il, ierr
416  INTEGER :: n_extrap                             ! Extrapolated points number
417  REAL, ALLOCATABLE :: ax(:), lon_rad(:), lon_ini(:), lev_dyn(:), yder(:)
418  REAL, ALLOCATABLE :: ay(:), lat_rad(:), lat_ini(:), var_tmp3d(:,:,:)
419  REAL, ALLOCATABLE, SAVE :: var_ana3d(:,:,:)
420!-------------------------------------------------------------------------------
421  iml=assert_eq(SIZE(lon_in),SIZE(lon_in2),SIZE(pls_in,1),SIZE(var3d,1),TRIM(modname)//" iml")
422  jml=assert_eq(SIZE(lat_in),              SIZE(pls_in,2),SIZE(var3d,2),TRIM(modname)//" jml")
423  lml=assert_eq(SIZE(pls_in,3),SIZE(var3d,3),TRIM(modname)//" lml"); jml2=SIZE(lat_in2)
424
425  WRITE(lunout, *)'Going into flinget to extract the 3D field.'
426  IF(.NOT.ALLOCATED(var_ana3d)) ALLOCATE(var_ana3d(iml_dyn, jml_dyn, llm_dyn))
427  CALL flinget(fid_dyn,var,iml_dyn,jml_dyn,llm_dyn,ttm_dyn,1,1,var_ana3d)
428
429!--- ANGLES IN DEGREES ARE CONVERTED INTO RADIANS
430  ALLOCATE(lon_ini(iml_dyn), lat_ini(jml_dyn))
431  lon_ini(:)=lon_dyn(:,1); IF(MAXVAL(lon_dyn)>pi) lon_ini=lon_ini*deg2rad
432  lat_ini(:)=lat_dyn(1,:); IF(MAXVAL(lat_dyn)>pi) lat_ini=lat_ini*deg2rad
433
434!--- FIELDS ARE PROCESSED TO BE ON STANDARD ANGULAR DOMAINS
435  ALLOCATE(lon_rad(iml_dyn), lat_rad(jml_dyn), lev_dyn(llm_dyn))
436  CALL conf_dat3d(var, lon_ini, lat_ini, levdyn_ini,                           &
437                       lon_rad, lat_rad, lev_dyn, var_ana3d, .TRUE.)
438  DEALLOCATE(lon_ini, lat_ini)
439
440!--- COMPUTE THE REQUIRED FIELDS USING ROUTINE grid_noro
441  ALLOCATE(var_tmp3d(iml,jml,llm_dyn))
442  DO il = 1,llm_dyn
443    CALL interp_startvar(var,il==1,lon_rad,lat_rad,var_ana3d(:,:,il),          &
444                          lon_in,lat_in,lon_in2,lat_in2,var_tmp3d(:,:,il))
445  END DO
446  DEALLOCATE(lon_rad, lat_rad)
447
448!--- VERTICAL INTERPOLATION FROM TOP OF ATMOSPHERE TO GROUND
449  ALLOCATE(ax(llm_dyn),ay(llm_dyn),yder(llm_dyn))
450  ax = lev_dyn(llm_dyn:1:-1)
451  skip = .FALSE.
452  n_extrap = 0
453  DO ij=1, jml
454    DO ii=1, iml-1
455      ay = var_tmp3d(ii, ij, llm_dyn:1:-1)
456      yder = pchsp_95(ax, ay, ibeg=2, iend=2, vc_beg=0., vc_end=0.)
457      CALL pchfe_95(ax, ay, yder, skip, pls_in(ii, ij, lml:1:-1),              &
458           var3d(ii, ij, lml:1:-1), ierr)
459      IF(ierr<0) CALL abort_gcm(TRIM(modname),'error in pchfe_95',1)
460      n_extrap = n_extrap + ierr
461    END DO
462  END DO
463  IF(n_extrap/=0) WRITE(lunout,*)TRIM(modname)//" pchfe_95: n_extrap=", n_extrap
464  var3d(iml, :, :) = var3d(1, :, :)
465
466  DO il=1, lml
467    CALL minmax(iml*jml, var3d(1, 1, il), chmin, chmax)
468    WRITE(lunout, *)' '//TRIM(var)//'  min max l ', il, chmin, chmax
469  END DO
470
471END SUBROUTINE start_inter_3d
472!
473!-------------------------------------------------------------------------------
474
475
476!-------------------------------------------------------------------------------
477!
478SUBROUTINE interp_startvar(nam,ibeg,lon,lat,vari,lon1,lat1,lon2,lat2,varo)
479!
480!-------------------------------------------------------------------------------
481  USE inter_barxy_m, ONLY: inter_barxy
482  IMPLICIT NONE
483!-------------------------------------------------------------------------------
484! Arguments:
485  CHARACTER(LEN=*), INTENT(IN)  :: nam
486  LOGICAL,          INTENT(IN)  :: ibeg
487  REAL,             INTENT(IN)  :: lon(:), lat(:)   ! dim (ii) (jj)
488  REAL,             INTENT(IN)  :: vari(:,:)        ! dim (ii,jj)
489  REAL,             INTENT(IN)  :: lon1(:), lat1(:) ! dim (i1) (j1)
490  REAL,             INTENT(IN)  :: lon2(:), lat2(:) ! dim (i1) (j2)
491  REAL,             INTENT(OUT) :: varo(:,:)        ! dim (i1) (j1)
492!-------------------------------------------------------------------------------
493! Local variables:
494  CHARACTER(LEN=256) :: modname="interp_startvar"
495  INTEGER            :: ii, jj, i1, j1, j2
496  REAL, ALLOCATABLE  :: vtmp(:,:)
497!-------------------------------------------------------------------------------
498  ii=assert_eq(SIZE(lon),            SIZE(vari,1),TRIM(modname)//" ii")
499  jj=assert_eq(SIZE(lat),            SIZE(vari,2),TRIM(modname)//" jj")
500  i1=assert_eq(SIZE(lon1),SIZE(lon2),SIZE(varo,1),TRIM(modname)//" i1")
501  j1=assert_eq(SIZE(lat1),           SIZE(varo,2),TRIM(modname)//" j1")
502  j2=SIZE(lat2)
503  ALLOCATE(vtmp(i1-1,j1))
504  IF(ibeg.AND.prt_level>1) THEN
505    WRITE(lunout,*)"---------------------------------------------------------"
506    WRITE(lunout,*)"$$$ Interpolation barycentrique pour "//TRIM(nam)//" $$$"
507    WRITE(lunout,*)"---------------------------------------------------------"
508  END IF
509  CALL inter_barxy(lon, lat(:jj-1), vari, lon2(:i1-1), lat2, vtmp)
510  CALL gr_int_dyn(vtmp, varo, i1-1, j1)
511
512END SUBROUTINE interp_startvar
513!
514!-------------------------------------------------------------------------------
515
516END MODULE etat0dyn
517!
518!*******************************************************************************
Note: See TracBrowser for help on using the repository browser.