source: LMDZ5/trunk/libf/dynphy_lonlat/phylmd/etat0dyn_netcdf.F90 @ 3094

Last change on this file since 3094 was 2941, checked in by acozic, 7 years ago

when we make a simulation with inca, we do not want to make a test on ozone for the ce0l

File size: 21.8 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 comconst_mod, ONLY: pi, cpp, kappa
39  USE comvert_mod, ONLY: ap, bp, preff, pressure_exner
40  USE temps_mod, ONLY: annee_ref, day_ref, itau_dyn, itau_phy, start_time
41 
42  IMPLICIT NONE
43
44  PRIVATE
45  PUBLIC :: etat0dyn_netcdf
46
47  include "iniprint.h"
48  include "dimensions.h"
49  include "paramet.h"
50  include "comgeom2.h"
51  include "comdissnew.h"
52  REAL, SAVE :: deg2rad
53  INTEGER,            SAVE      :: iml_dyn, jml_dyn, llm_dyn, ttm_dyn, fid_dyn
54  REAL, ALLOCATABLE,  SAVE      :: lon_dyn(:,:), lat_dyn(:,:), levdyn_ini(:)
55  CHARACTER(LEN=120), PARAMETER :: dynfname='ECDYN.nc'
56
57CONTAINS
58
59!-------------------------------------------------------------------------------
60!
61SUBROUTINE etat0dyn_netcdf(masque, phis)
62!
63!-------------------------------------------------------------------------------
64! Purpose: Create dynamical initial states.
65!-------------------------------------------------------------------------------
66! Notes:  1) This routine is designed to work for Earth
67!         2) If masque(:,:)/=-99999., masque and phis are already known.
68!         Otherwise: compute it.
69!-------------------------------------------------------------------------------
70  USE control_mod
71  USE regr_lat_time_coefoz_m, ONLY: regr_lat_time_coefoz
72  USE regr_pr_o3_m,   ONLY: regr_pr_o3
73  USE press_coefoz_m, ONLY: press_coefoz
74  USE exner_hyb_m,    ONLY: exner_hyb
75  USE exner_milieu_m, ONLY: exner_milieu
76  USE infotrac,       ONLY: nqtot, tname
77  USE filtreg_mod
78  IMPLICIT NONE
79!-------------------------------------------------------------------------------
80! Arguments:
81  REAL,    INTENT(INOUT) :: masque(iip1,jjp1)   !--- Land-ocean mask
82  REAL,    INTENT(INOUT) :: phis  (iip1,jjp1)   !--- Ground geopotential
83!-------------------------------------------------------------------------------
84! Local variables:
85  CHARACTER(LEN=256) :: modname, fmt
86  INTEGER            :: i, j, l, ji, itau, iday
87  REAL               :: xpn, xps, time, phystep
88  REAL, DIMENSION(iip1,jjp1)       :: psol
89  REAL, DIMENSION(iip1,jjp1,llm+1) :: p3d
90  REAL, DIMENSION(iip1,jjp1,llm)   :: uvent, t3d, tpot, qsat, qd
91  REAL, DIMENSION(iip1,jjp1,llm)   :: pk, pls, y, masse
92  REAL, DIMENSION(iip1,jjm ,llm)   :: vvent
93  REAL, DIMENSION(ip1jm    ,llm)   :: pbarv
94  REAL, DIMENSION(ip1jmp1  ,llm)   :: pbaru, phi, w
95  REAL, DIMENSION(ip1jmp1)         :: pks
96  REAL, DIMENSION(iim)             :: xppn, xpps
97  REAL, ALLOCATABLE                :: q3d(:,:,:,:)
98!-------------------------------------------------------------------------------
99  modname='etat0dyn_netcdf'
100
101  deg2rad = pi/180.0
102
103! Compute psol AND tsol, knowing phis.
104!*******************************************************************************
105  CALL start_init_dyn(rlonv, rlatu, rlonu, rlatv, phis, psol)
106
107! Mid-levels pressure computation
108!*******************************************************************************
109  CALL pression(ip1jmp1, ap, bp, psol, p3d)             !--- Update p3d
110  IF(pressure_exner) THEN                               !--- Update pk, pks
111    CALL exner_hyb   (ip1jmp1,psol,p3d,pks,pk)
112  ELSE
113    CALL exner_milieu(ip1jmp1,psol,p3d,pks,pk)
114  END IF
115  pls(:,:,:)=preff*(pk(:,:,:)/cpp)**(1./kappa)          !--- Update pls
116
117! Update uvent, vvent, t3d and tpot
118!*******************************************************************************
119  uvent(:,:,:) = 0.0 ; vvent(:,:,:) = 0.0 ; t3d (:,:,:) = 0.0
120  CALL startget_dyn3d('u'   ,rlonu,rlatu,pls,y ,uvent,rlonv,rlatv)
121  CALL startget_dyn3d('v'   ,rlonv,rlatv,pls(:,:jjm,:),y(:,:jjm,:),vvent,      &
122 &                           rlonu,rlatu(:jjm))
123  CALL startget_dyn3d('t'   ,rlonv,rlatu,pls,y ,t3d ,rlonu,rlatv)
124  tpot(:,:,:)=t3d(:,:,:)
125  CALL startget_dyn3d('tpot',rlonv,rlatu,pls,pk,tpot,rlonu,rlatv)
126
127  WRITE(lunout,*) 'T3D min,max:',MINVAL(t3d(:,:,:)),MAXVAL(t3d(:,:,:))
128  WRITE(lunout,*) 'PLS min,max:',MINVAL(pls(:,:,:)),MAXVAL(pls(:,:,:))
129
130! Humidity at saturation computation
131!*******************************************************************************
132  WRITE(lunout,*) 'avant q_sat'
133  CALL q_sat(llm*jjp1*iip1, t3d, pls, qsat)
134  WRITE(lunout,*) 'apres q_sat'
135  WRITE(lunout,*) 'QSAT min,max:',MINVAL(qsat(:,:,:)),MAXVAL(qsat(:,:,:))
136!  WRITE(lunout,*) 'QSAT :',qsat(10,20,:)
137  qd (:,:,:) = 0.0
138  CALL startget_dyn3d('q',rlonv,rlatu,pls,qsat,qd,rlonu,rlatv)
139  ALLOCATE(q3d(iip1,jjp1,llm,nqtot)); q3d(:,:,:,:)=0.0 ; q3d(:,:,:,1)=qd(:,:,:)
140  CALL flinclo(fid_dyn)
141
142! Parameterization of ozone chemistry:
143!*******************************************************************************
144! Look for ozone tracer:
145#ifndef INCA
146  DO i=1,nqtot; IF(ANY(["O3","o3"]==tname(i))) EXIT; END DO
147  IF(i/=nqtot+1) THEN
148    CALL regr_lat_time_coefoz
149    CALL press_coefoz
150    CALL regr_pr_o3(p3d, q3d(:,:,:,i))
151    q3d(:,:,:,i)=q3d(:,:,:,i)*48./ 29.                  !--- Mole->mass fraction         
152  END IF
153#endif
154  q3d(iip1,:,:,:)=q3d(1,:,:,:)
155
156! Writing
157!*******************************************************************************
158  CALL inidissip(lstardis, nitergdiv, nitergrot, niterh, tetagdiv, tetagrot,   &
159       tetatemp, vert_prof_dissip)
160  WRITE(lunout,*)'sortie inidissip'
161  itau=0
162  itau_dyn=0
163  itau_phy=0
164  iday=dayref+itau/day_step
165  time=FLOAT(itau-(iday-dayref)*day_step)/day_step
166  IF(time>1.) THEN
167   time=time-1
168   iday=iday+1
169  END IF
170  day_ref=dayref
171  annee_ref=anneeref
172  CALL geopot( ip1jmp1, tpot, pk, pks, phis, phi )
173  WRITE(lunout,*)'sortie geopot'
174  CALL caldyn0( itau, uvent, vvent, tpot, psol, masse, pk, phis,               &
175                phi,  w, pbaru, pbarv, time+iday-dayref)
176  WRITE(lunout,*)'sortie caldyn0'
177  start_time = 0.
178#ifdef CPP_PARA
179  CALL dynredem0_loc( "start.nc", dayref, phis)
180#else
181  CALL dynredem0( "start.nc", dayref, phis)
182#endif
183  WRITE(lunout,*)'sortie dynredem0'
184#ifdef CPP_PARA
185  CALL dynredem1_loc( "start.nc", 0.0, vvent, uvent, tpot, q3d, masse, psol)
186#else
187  CALL dynredem1( "start.nc", 0.0, vvent, uvent, tpot, q3d, masse, psol)
188#endif
189  WRITE(lunout,*)'sortie dynredem1'
190  CALL histclo()
191
192END SUBROUTINE etat0dyn_netcdf
193!
194!-------------------------------------------------------------------------------
195
196
197!-------------------------------------------------------------------------------
198!
199SUBROUTINE startget_dyn3d(var, lon_in,  lat_in,  pls,  workvar,&
200                        champ, lon_in2, lat_in2)
201!-------------------------------------------------------------------------------
202  IMPLICIT NONE
203!===============================================================================
204! Purpose: Compute some quantities (u,v,t,q,tpot) using variables U,V,TEMP and R
205!     (3D fields) of file dynfname.
206!-------------------------------------------------------------------------------
207! Note: An input auxilliary field "workvar" has to be specified in two cases:
208!     * for "q":    the saturated humidity.
209!     * for "tpot": the Exner function.
210!===============================================================================
211! Arguments:
212  CHARACTER(LEN=*), INTENT(IN)    :: var
213  REAL,             INTENT(IN)    :: lon_in(:)        ! dim (iml)
214  REAL,             INTENT(IN)    :: lat_in(:)        ! dim (jml)
215  REAL,             INTENT(IN)    :: pls    (:, :, :) ! dim (iml, jml, lml)
216  REAL,             INTENT(IN)    :: workvar(:, :, :) ! dim (iml, jml, lml)
217  REAL,             INTENT(INOUT) :: champ  (:, :, :) ! dim (iml, jml, lml)
218  REAL,             INTENT(IN)    :: lon_in2(:)       ! dim (iml)
219  REAL,             INTENT(IN)    :: lat_in2(:)       ! dim (jml2)
220!-------------------------------------------------------------------------------
221! Local variables:
222  CHARACTER(LEN=10)  :: vname
223  CHARACTER(LEN=256) :: msg, modname="startget_dyn3d"
224  INTEGER            :: iml, jml, jml2, lml, il
225  REAL               :: xppn, xpps
226!-------------------------------------------------------------------------------
227  iml=assert_eq([SIZE(lon_in),SIZE(pls,1),SIZE(workvar,1),SIZE(champ,1),       &
228     &                                    SIZE(lon_in2)], TRIM(modname)//" iml")
229  jml=assert_eq( SIZE(lat_in),SIZE(pls,2),SIZE(workvar,2),SIZE(champ,2),       &
230     &                                                    TRIM(modname)//" jml")
231  lml=assert_eq(              SIZE(pls,3),SIZE(workvar,3),SIZE(champ,3),       &
232     &                                                    TRIM(modname)//" lml")
233  jml2=SIZE(lat_in2)
234
235!--- CHECK IF THE FIELD IS KNOWN
236   SELECT CASE(var)
237    CASE('u');    vname='U'
238    CASE('v');    vname='V'
239    CASE('t');    vname='TEMP'
240    CASE('q');    vname='R';    msg='humidity as the saturated humidity'
241    CASE('tpot'); msg='potential temperature as the Exner function'
242    CASE DEFAULT; msg='No rule to extract variable '//TRIM(var)
243      CALL abort_gcm(modname,TRIM(msg)//' from any data set',1)
244  END SELECT
245
246!--- CHECK IF SOMETHING IS MISSING
247  IF((var=='tpot'.OR.var=='q').AND.MINVAL(workvar)==MAXVAL(workvar)) THEN
248    msg='Could not compute '//TRIM(msg)//' is missing or constant.'
249    CALL abort_gcm(modname,TRIM(msg),1)
250  END IF
251
252!--- INTERPOLATE 3D FIELD IF NEEDED
253  IF(var/='tpot') CALL start_inter_3d(TRIM(vname),lon_in,lat_in,lon_in2,      &
254                                                  lat_in2,pls,champ)
255
256!--- COMPUTE THE REQUIRED FILED
257  SELECT CASE(var)
258    CASE('u'); DO il=1,lml; champ(:,:,il)=champ(:,:,il)*cu(:,1:jml); END DO
259      champ(iml,:,:)=champ(1,:,:)                   !--- Eastward wind
260
261    CASE('v'); DO il=1,lml; champ(:,:,il)=champ(:,:,il)*cv(:,1:jml); END DO
262      champ(iml,:,:)=champ(1,:,:)                   !--- Northward wind
263
264    CASE('tpot','q')
265      IF(var=='tpot') THEN; champ=champ*cpp/workvar !--- Potential temperature
266      ELSE;                 champ=champ*.01*workvar !--- Relative humidity
267        WHERE(champ<0.) champ=1.0E-10
268      END IF
269      DO il=1,lml
270        xppn = SUM(aire(:,1  )*champ(:,1  ,il))/apoln
271        xpps = SUM(aire(:,jml)*champ(:,jml,il))/apols
272        champ(:,1  ,il) = xppn
273        champ(:,jml,il) = xpps
274      END DO
275  END SELECT
276
277END SUBROUTINE startget_dyn3d
278!
279!-------------------------------------------------------------------------------
280
281
282!-------------------------------------------------------------------------------
283!
284SUBROUTINE start_init_dyn(lon_in,lat_in,lon_in2,lat_in2,zs,psol)
285!
286!-------------------------------------------------------------------------------
287  IMPLICIT NONE
288!===============================================================================
289! Purpose:   Compute psol, knowing phis.
290!===============================================================================
291! Arguments:
292  REAL,    INTENT(IN)  :: lon_in (:),  lat_in (:)    ! dim (iml) (jml)
293  REAL,    INTENT(IN)  :: lon_in2(:),  lat_in2(:)    ! dim (iml) (jml2)
294  REAL,    INTENT(IN)  :: zs  (:,:)                  ! dim (iml,jml)
295  REAL,    INTENT(OUT) :: psol(:,:)                  ! dim (iml,jml)
296!-------------------------------------------------------------------------------
297! Local variables:
298  CHARACTER(LEN=256) :: modname='start_init_dyn'
299  REAL               :: date, dt
300  INTEGER            :: iml, jml, jml2, itau(1)
301  REAL, ALLOCATABLE  :: lon_rad(:), lon_ini(:), var_ana(:,:)
302  REAL, ALLOCATABLE  :: lat_rad(:), lat_ini(:)
303  REAL, ALLOCATABLE  :: z(:,:), ps(:,:), ts(:,:)
304!-------------------------------------------------------------------------------
305  iml=assert_eq(SIZE(lon_in),SIZE(zs,1),SIZE(psol,1),SIZE(lon_in2),            &
306      &                                              TRIM(modname)//" iml")
307  jml=assert_eq(SIZE(lat_in),SIZE(zs,2),SIZE(psol,2),TRIM(modname)//" jml")
308  jml2=SIZE(lat_in2)
309
310  WRITE(lunout,*) 'Opening the surface analysis'
311  CALL flininfo(dynfname, iml_dyn, jml_dyn, llm_dyn, ttm_dyn, fid_dyn)
312  WRITE(lunout,*) 'Values read: ', iml_dyn, jml_dyn, llm_dyn, ttm_dyn
313
314  ALLOCATE(lon_dyn(iml_dyn,jml_dyn), lat_dyn(iml_dyn,jml_dyn))
315  ALLOCATE(levdyn_ini(llm_dyn))
316  CALL flinopen(dynfname, .FALSE., iml_dyn, jml_dyn, llm_dyn,                  &
317                lon_dyn,lat_dyn,levdyn_ini,ttm_dyn,itau,date,dt,fid_dyn)
318
319!--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS
320  ALLOCATE(lon_ini(iml_dyn),lat_ini(jml_dyn))
321  lon_ini(:)=lon_dyn(:,1); IF(MAXVAL(lon_dyn)>pi) lon_ini=lon_ini*deg2rad
322  lat_ini(:)=lat_dyn(1,:); IF(MAXVAL(lat_dyn)>pi) lat_ini=lat_ini*deg2rad
323
324  ALLOCATE(var_ana(iml_dyn,jml_dyn),lon_rad(iml_dyn),lat_rad(jml_dyn))
325  CALL get_var_dyn('Z',z)                        !--- SURFACE GEOPOTENTIAL
326  CALL get_var_dyn('SP',ps)                      !--- SURFACE PRESSURE
327  CALL get_var_dyn('ST',ts)                      !--- SURFACE TEMPERATURE
328!  CALL flinclo(fid_dyn)
329  DEALLOCATE(var_ana,lon_rad,lat_rad,lon_ini,lat_ini)
330
331!--- PSOL IS COMPUTED IN PASCALS
332  psol(:iml-1,:) = ps(:iml-1,:)*(1.0+(z(:iml-1,:)-zs(:iml-1,:))/287.0          &
333     &            /ts(:iml-1,:))
334  psol(iml,:)=psol(1,:)
335  DEALLOCATE(z,ps,ts)
336  psol(:,1  )=SUM(aire(1:iml-1,1  )*psol(1:iml-1,1  ))/apoln  !--- NORTH POLE
337  psol(:,jml)=SUM(aire(1:iml-1,jml)*psol(1:iml-1,jml))/apols  !--- SOUTH POLE
338
339CONTAINS
340
341!-------------------------------------------------------------------------------
342!
343SUBROUTINE get_var_dyn(title,field)
344!
345!-------------------------------------------------------------------------------
346  USE conf_dat_m, ONLY: conf_dat2d
347  IMPLICIT NONE
348!-------------------------------------------------------------------------------
349! Arguments:
350  CHARACTER(LEN=*),  INTENT(IN)    :: title
351  REAL, ALLOCATABLE, INTENT(INOUT) :: field(:,:)
352!-------------------------------------------------------------------------------
353! Local variables:
354  CHARACTER(LEN=256) :: msg
355  INTEGER :: tllm
356!-------------------------------------------------------------------------------
357  SELECT CASE(title)
358    CASE('Z');     tllm=0;       msg='geopotential'
359    CASE('SP');    tllm=0;       msg='surface pressure'
360    CASE('ST');    tllm=llm_dyn; msg='temperature'
361  END SELECT
362  IF(.NOT.ALLOCATED(field)) THEN
363    ALLOCATE(field(iml,jml))
364    CALL flinget(fid_dyn, title, iml_dyn,jml_dyn, tllm, ttm_dyn, 1, 1, var_ana)
365    CALL conf_dat2d(title, lon_ini, lat_ini, lon_rad, lat_rad, var_ana, .TRUE.)
366    CALL interp_startvar(title, .TRUE., lon_rad,lat_rad, var_ana,              &
367                                        lon_in, lat_in, lon_in2, lat_in2, field)
368  ELSE IF(SIZE(field)/=SIZE(z)) THEN
369    msg='The '//TRIM(msg)//' field we have does not have the right size'
370    CALL abort_gcm(TRIM(modname),msg,1)
371  END IF
372
373END SUBROUTINE get_var_dyn
374!
375!-------------------------------------------------------------------------------
376
377END SUBROUTINE start_init_dyn
378!
379!-------------------------------------------------------------------------------
380
381
382!-------------------------------------------------------------------------------
383!
384SUBROUTINE start_inter_3d(var,lon_in,lat_in,lon_in2,lat_in2,pls_in,var3d)
385!
386!-------------------------------------------------------------------------------
387  USE conf_dat_m, ONLY: conf_dat3d
388  USE pchsp_95_m, ONLY: pchsp_95
389  USE pchfe_95_m, ONLY: pchfe_95
390  IMPLICIT NONE
391!-------------------------------------------------------------------------------
392! Arguments:
393  CHARACTER(LEN=*), INTENT(IN) :: var
394  REAL,    INTENT(IN)  :: lon_in(:),  lat_in(:)   ! dim (iml) (jml)
395  REAL,    INTENT(IN)  :: lon_in2(:), lat_in2(:)  ! dim (iml) (jml2)
396  REAL,    INTENT(IN)  :: pls_in(:,:,:)           ! dim (iml,jml,lml)
397  REAL,    INTENT(OUT) :: var3d (:,:,:)           ! dim (iml,jml,lml)
398!-------------------------------------------------------------------------------
399! Local variables:
400  CHARACTER(LEN=256) :: modname='start_inter_3d'
401  LOGICAL :: skip
402  REAL    :: chmin, chmax
403  INTEGER :: iml, jml, lml, jml2, ii, ij, il, ierr
404  INTEGER :: n_extrap                             ! Extrapolated points number
405  REAL, ALLOCATABLE :: ax(:), lon_rad(:), lon_ini(:), lev_dyn(:), yder(:)
406  REAL, ALLOCATABLE :: ay(:), lat_rad(:), lat_ini(:), var_tmp3d(:,:,:)
407  REAL, ALLOCATABLE, SAVE :: var_ana3d(:,:,:)
408!-------------------------------------------------------------------------------
409  iml=assert_eq(SIZE(lon_in),SIZE(lon_in2),SIZE(pls_in,1),SIZE(var3d,1),TRIM(modname)//" iml")
410  jml=assert_eq(SIZE(lat_in),              SIZE(pls_in,2),SIZE(var3d,2),TRIM(modname)//" jml")
411  lml=assert_eq(SIZE(pls_in,3),SIZE(var3d,3),TRIM(modname)//" lml"); jml2=SIZE(lat_in2)
412
413  WRITE(lunout, *)'Going into flinget to extract the 3D field.'
414  IF(.NOT.ALLOCATED(var_ana3d)) ALLOCATE(var_ana3d(iml_dyn, jml_dyn, llm_dyn))
415  CALL flinget(fid_dyn,var,iml_dyn,jml_dyn,llm_dyn,ttm_dyn,1,1,var_ana3d)
416
417!--- ANGLES IN DEGREES ARE CONVERTED INTO RADIANS
418  ALLOCATE(lon_ini(iml_dyn), lat_ini(jml_dyn))
419  lon_ini(:)=lon_dyn(:,1); IF(MAXVAL(lon_dyn)>pi) lon_ini=lon_ini*deg2rad
420  lat_ini(:)=lat_dyn(1,:); IF(MAXVAL(lat_dyn)>pi) lat_ini=lat_ini*deg2rad
421
422!--- FIELDS ARE PROCESSED TO BE ON STANDARD ANGULAR DOMAINS
423  ALLOCATE(lon_rad(iml_dyn), lat_rad(jml_dyn), lev_dyn(llm_dyn))
424  CALL conf_dat3d(var, lon_ini, lat_ini, levdyn_ini,                           &
425                       lon_rad, lat_rad, lev_dyn, var_ana3d, .TRUE.)
426  DEALLOCATE(lon_ini, lat_ini)
427
428!--- COMPUTE THE REQUIRED FIELDS USING ROUTINE grid_noro
429  ALLOCATE(var_tmp3d(iml,jml,llm_dyn))
430  DO il = 1,llm_dyn
431    CALL interp_startvar(var,il==1,lon_rad,lat_rad,var_ana3d(:,:,il),          &
432                          lon_in,lat_in,lon_in2,lat_in2,var_tmp3d(:,:,il))
433  END DO
434  DEALLOCATE(lon_rad, lat_rad)
435
436!--- VERTICAL INTERPOLATION FROM TOP OF ATMOSPHERE TO GROUND
437  ALLOCATE(ax(llm_dyn),ay(llm_dyn),yder(llm_dyn))
438  ax = lev_dyn(llm_dyn:1:-1)
439  skip = .FALSE.
440  n_extrap = 0
441  DO ij=1, jml
442    DO ii=1, iml-1
443      ay = var_tmp3d(ii, ij, llm_dyn:1:-1)
444      yder = pchsp_95(ax, ay, ibeg=2, iend=2, vc_beg=0., vc_end=0.)
445      CALL pchfe_95(ax, ay, yder, skip, pls_in(ii, ij, lml:1:-1),              &
446           var3d(ii, ij, lml:1:-1), ierr)
447      IF(ierr<0) CALL abort_gcm(TRIM(modname),'error in pchfe_95',1)
448      n_extrap = n_extrap + ierr
449    END DO
450  END DO
451  IF(n_extrap/=0) WRITE(lunout,*)TRIM(modname)//" pchfe_95: n_extrap=", n_extrap
452  var3d(iml, :, :) = var3d(1, :, :)
453
454  DO il=1, lml
455    CALL minmax(iml*jml, var3d(1, 1, il), chmin, chmax)
456    WRITE(lunout, *)' '//TRIM(var)//'  min max l ', il, chmin, chmax
457  END DO
458
459END SUBROUTINE start_inter_3d
460!
461!-------------------------------------------------------------------------------
462
463
464!-------------------------------------------------------------------------------
465!
466SUBROUTINE interp_startvar(nam,ibeg,lon,lat,vari,lon1,lat1,lon2,lat2,varo)
467!
468!-------------------------------------------------------------------------------
469  USE inter_barxy_m, ONLY: inter_barxy
470  IMPLICIT NONE
471!-------------------------------------------------------------------------------
472! Arguments:
473  CHARACTER(LEN=*), INTENT(IN)  :: nam
474  LOGICAL,          INTENT(IN)  :: ibeg
475  REAL,             INTENT(IN)  :: lon(:), lat(:)   ! dim (ii) (jj)
476  REAL,             INTENT(IN)  :: vari(:,:)        ! dim (ii,jj)
477  REAL,             INTENT(IN)  :: lon1(:), lat1(:) ! dim (i1) (j1)
478  REAL,             INTENT(IN)  :: lon2(:), lat2(:) ! dim (i1) (j2)
479  REAL,             INTENT(OUT) :: varo(:,:)        ! dim (i1) (j1)
480!-------------------------------------------------------------------------------
481! Local variables:
482  CHARACTER(LEN=256) :: modname="interp_startvar"
483  INTEGER            :: ii, jj, i1, j1, j2
484  REAL, ALLOCATABLE  :: vtmp(:,:)
485!-------------------------------------------------------------------------------
486  ii=assert_eq(SIZE(lon),            SIZE(vari,1),TRIM(modname)//" ii")
487  jj=assert_eq(SIZE(lat),            SIZE(vari,2),TRIM(modname)//" jj")
488  i1=assert_eq(SIZE(lon1),SIZE(lon2),SIZE(varo,1),TRIM(modname)//" i1")
489  j1=assert_eq(SIZE(lat1),           SIZE(varo,2),TRIM(modname)//" j1")
490  j2=SIZE(lat2)
491  ALLOCATE(vtmp(i1-1,j1))
492  IF(ibeg.AND.prt_level>1) THEN
493    WRITE(lunout,*)"---------------------------------------------------------"
494    WRITE(lunout,*)"$$$ Interpolation barycentrique pour "//TRIM(nam)//" $$$"
495    WRITE(lunout,*)"---------------------------------------------------------"
496  END IF
497  CALL inter_barxy(lon, lat(:jj-1), vari, lon2(:i1-1), lat2, vtmp)
498  CALL gr_int_dyn(vtmp, varo, i1-1, j1)
499
500END SUBROUTINE interp_startvar
501!
502!-------------------------------------------------------------------------------
503
504END MODULE etat0dyn
505!
506!*******************************************************************************
Note: See TracBrowser for help on using the repository browser.