source: LMDZ6/branches/Amaury_dev/libf/dynphy_lonlat/phylmd/etat0dyn_netcdf.F90 @ 5117

Last change on this file since 5117 was 5117, checked in by abarral, 4 months ago

rename modules properly lmdz_*
move some unused files to obsolete/
(lint) uppercase fortran keywords

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