source: LMDZ5/trunk/libf/dynlonlat_phylonlat/phylmd/etat0phys_netcdf.F90 @ 2396

Last change on this file since 2396 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: 18.9 KB
Line 
1MODULE etat0phys
2!
3!*******************************************************************************
4! Purpose: Create physical 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!    *  etat0phys_netcdf routine can access to NetCDF data through subroutines:
12!         "start_init_phys" for variables contained in file "ECPHY.nc":
13!            'ST'     : Surface temperature
14!            'CDSW'   : Soil moisture
15!         "start_init_orog" for variables contained in file "Relief.nc":
16!            'RELIEF' : High resolution orography
17!
18!    * The land mask and corresponding weights can be:
19!      1) computed using the ocean mask from the ocean model (to ensure ocean
20!         fractions are the same for atmosphere and ocean) for coupled runs.
21!         File name: "o2a.nc"  ;  variable name: "OceMask"
22!      2) computed from topography file "Relief.nc" for forced runs.
23!
24!    * Allowed values for read_climoz flag are 0, 1 and 2:
25!      0: do not read an ozone climatology
26!      1: read a single ozone climatology that will be used day and night
27!      2: read two ozone climatologies, the average day and night climatology
28!         and the daylight climatology
29!-------------------------------------------------------------------------------
30!    * There is a big mess with the longitude size. Should it be iml or iml+1 ?
31!  I have chosen to use the iml+1 as an argument to this routine and we declare
32!  internaly smaller fields when needed. This needs to be cleared once and for
33!  all in LMDZ. A convention is required.
34!-------------------------------------------------------------------------------
35
36  USE ioipsl,             ONLY: flininfo, flinopen, flinget, flinclo
37  USE assert_eq_m,        ONLY: assert_eq
38  USE dimphy,             ONLY: klon
39  USE conf_dat_m,         ONLY: conf_dat2d
40  USE phys_state_var_mod, ONLY: zmea, zstd, zsig, zgam, zthe, zpic, zval, z0m, &
41    rlon, solsw, radsol, t_ancien, wake_deltat, wake_s,  rain_fall, qsol, z0h, &
42    rlat, sollw, rugoro, q_ancien, wake_deltaq, wake_pe, snow_fall, ratqs,w01, &
43    sig1, ftsol, clwcon, fm_therm, wake_Cstar,  pctsrf,  entr_therm,radpas, f0,&
44    zmax0,fevap, rnebcon,falb_dir, wake_fip,    agesno,  detr_therm, pbl_tke,  &
45    phys_state_var_init
46
47  PRIVATE
48  PUBLIC :: etat0phys_netcdf
49
50  include "iniprint.h"
51  include "dimensions.h"
52  include "paramet.h"
53  include "comgeom2.h"
54  include "comconst.h"
55  include "dimsoil.h"
56  include "temps.h"
57  include "clesphys.h"
58  REAL, SAVE :: deg2rad
59  REAL, SAVE, ALLOCATABLE :: tsol(:)
60  INTEGER,            SAVE      :: iml_phys, jml_phys, llm_phys, ttm_phys, fid_phys
61  REAL, ALLOCATABLE,  SAVE      :: lon_phys(:,:), lat_phys(:,:), levphys_ini(:)
62  CHARACTER(LEN=256), PARAMETER :: orofname="Relief.nc", orogvar="RELIEF"
63  CHARACTER(LEN=256), PARAMETER :: phyfname="ECPHY.nc",  psrfvar="SP"
64  CHARACTER(LEN=256), PARAMETER :: qsolvar="CDSW",       tsrfvar="ST"
65
66
67CONTAINS
68
69
70!-------------------------------------------------------------------------------
71!
72SUBROUTINE etat0phys_netcdf(masque, phis)
73!
74!-------------------------------------------------------------------------------
75! Purpose: Creates initial states
76!-------------------------------------------------------------------------------
77! Notes:  1) This routine is designed to work for Earth
78!         2) If masque(:,:)/=-99999., masque and phis are already known.
79!         Otherwise: compute it.
80!-------------------------------------------------------------------------------
81  USE control_mod
82  USE fonte_neige_mod
83  USE pbl_surface_mod
84  USE regr_lat_time_climoz_m, ONLY: regr_lat_time_climoz
85  USE indice_sol_mod
86  USE conf_phys_m, ONLY: conf_phys
87  USE init_ssrf_m, ONLY: start_init_subsurf
88  IMPLICIT NONE
89!-------------------------------------------------------------------------------
90! Arguments:
91  REAL,    INTENT(INOUT) :: masque(:,:) !--- Land mask           dim(iip1,jjp1)
92  REAL,    INTENT(INOUT) :: phis  (:,:) !--- Ground geopotential dim(iip1,jjp1)
93!-------------------------------------------------------------------------------
94! Local variables:
95  CHARACTER(LEN=256) :: modname="etat0phys_netcdf", fmt
96  INTEGER            :: i, j, l, ji, iml, jml
97  LOGICAL            :: read_mask
98  REAL               :: phystep, dummy
99  REAL, DIMENSION(SIZE(masque,1),SIZE(masque,2)) :: masque_tmp
100  REAL, DIMENSION(klon)               :: sn, rugmer, run_off_lic_0, fder
101  REAL, DIMENSION(klon,nbsrf)         :: qsolsrf, snsrf
102  REAL, DIMENSION(klon,nsoilmx,nbsrf) :: tsoil
103
104!--- Arguments for conf_phys
105  LOGICAL :: ok_journe, ok_mensuel, ok_instan, ok_hf, ok_LES, callstats
106  REAL    :: solarlong0, seuil_inversion, fact_cldcon, facttemps
107  LOGICAL :: ok_newmicro
108  INTEGER :: iflag_radia, iflag_cldcon, iflag_ratqs
109  REAL    :: ratqsbas, ratqshaut, tau_ratqs
110  LOGICAL :: ok_ade, ok_aie, ok_cdnc, aerosol_couple
111  INTEGER :: flag_aerosol
112  LOGICAL :: flag_aerosol_strat
113  LOGICAL :: new_aod
114  REAL    :: bl95_b0, bl95_b1
115  INTEGER :: read_climoz                        !--- Read ozone climatology
116  REAL    :: alp_offset
117
118  deg2rad= pi/180.0
119  iml=assert_eq(SIZE(masque,1),SIZE(phis,1),TRIM(modname)//" iml")
120  jml=assert_eq(SIZE(masque,2),SIZE(phis,2),TRIM(modname)//" jml")
121
122! Physics configuration
123!*******************************************************************************
124  CALL conf_phys(  ok_journe, ok_mensuel, ok_instan, ok_hf, ok_LES,     &
125                   callstats,                                           &
126                   solarlong0,seuil_inversion,                          &
127                   fact_cldcon, facttemps,ok_newmicro,iflag_radia,      &
128                   iflag_cldcon,                                        &
129                   iflag_ratqs,ratqsbas,ratqshaut,tau_ratqs,            &
130                   ok_ade, ok_aie, ok_cdnc, aerosol_couple,             &
131                   flag_aerosol, flag_aerosol_strat, new_aod,           &
132                   bl95_b0, bl95_b1,                                    &
133                   read_climoz,                                         &
134                   alp_offset)
135  CALL phys_state_var_init(read_climoz)
136
137!--- Initial atmospheric CO2 conc. from .def file
138  co2_ppm0 = co2_ppm
139
140! Compute ground geopotential, sub-cells quantities and possibly the mask.
141!*******************************************************************************
142  read_mask=ANY(masque/=-99999.); masque_tmp=masque
143  CALL start_init_orog(rlonv, rlatu, phis, masque_tmp)
144  WRITE(fmt,"(i4,'i1)')")iml ; fmt='('//ADJUSTL(fmt)
145  IF(.NOT.read_mask) THEN                       !--- Keep mask form orography
146    masque=masque_tmp
147    IF(prt_level>=1) THEN
148      WRITE(lunout,*)'BUILT MASK :'
149      WRITE(lunout,fmt) NINT(masque)
150    END IF
151    WHERE(   masque(:,:)<EPSFRA) masque(:,:)=0.
152    WHERE(1.-masque(:,:)<EPSFRA) masque(:,:)=1.
153  END IF
154  CALL gr_dyn_fi(1,iml,jml,klon,masque,zmasq) !--- Land mask to physical grid
155
156! Compute tsol and qsol on physical grid, knowing phis on 2D grid.
157!*******************************************************************************
158  CALL start_init_phys(rlonu, rlatv, phis)
159
160! Some initializations.
161!*******************************************************************************
162  sn    (:) = 0.0                               !--- Snow
163  radsol(:) = 0.0                               !--- Net radiation at ground
164  rugmer(:) = 0.001                             !--- Ocean rugosity
165  IF(read_climoz>=1) &                          !--- Ozone climatology
166    CALL regr_lat_time_climoz(read_climoz)
167
168! Sub-surfaces initialization.
169!*******************************************************************************
170  CALL start_init_subsurf(read_mask)
171
172! Write physical initial state
173!*******************************************************************************
174  WRITE(lunout,*)'phystep ',dtvr,iphysiq,nbapp_rad
175  phystep = dtvr * FLOAT(iphysiq)
176  radpas  = NINT (86400./phystep/ FLOAT(nbapp_rad) )
177  WRITE(lunout,*)'phystep =', phystep, radpas
178
179! Init: ftsol, snsrf, qsolsrf, tsoil, rain_fall, snow_fall, solsw, sollw, z0
180!*******************************************************************************
181  DO i=1,nbsrf; ftsol(:,i) = tsol; END DO
182  DO i=1,nbsrf; snsrf(:,i) = sn;   END DO
183  falb_dir(:,is_ter,:) = 0.08
184  falb_dir(:,is_lic,:) = 0.6
185  falb_dir(:,is_oce,:) = 0.5
186  falb_dir(:,is_sic,:) = 0.6
187  fevap(:,:) = 0.
188  DO i=1,nbsrf; qsolsrf(:,i)=150.; END DO
189  DO i=1,nbsrf; DO j=1,nsoilmx; tsoil(:,j,i) = tsol; END DO; END DO
190  rain_fall  = 0.
191  snow_fall  = 0.
192  solsw      = 165.
193  sollw      = -53.
194  t_ancien   = 273.15
195  q_ancien   = 0.
196  agesno     = 0.
197
198  z0m(:,is_oce) = rugmer(:)
199  z0m(:,is_ter) = MAX(1.0e-05,zstd(:)*zsig(:)/2.0)
200  z0m(:,is_lic) = MAX(1.0e-05,zstd(:)*zsig(:)/2.0)
201  z0m(:,is_sic) = 0.001
202  z0h(:,:)=z0m(:,:)
203
204  fder    = 0.0
205  clwcon  = 0.0
206  rnebcon = 0.0
207  ratqs   = 0.0
208  run_off_lic_0 = 0.0
209  rugoro  = 0.0
210
211! Before phyredem calling, surface modules and values to be saved in startphy.nc
212! are initialized
213!*******************************************************************************
214  dummy            = 1.0
215  pbl_tke(:,:,:)   = 1.e-8
216  zmax0(:)         = 40.
217  f0(:)            = 1.e-5
218  sig1(:,:)        = 0.
219  w01(:,:)         = 0.
220  wake_deltat(:,:) = 0.
221  wake_deltaq(:,:) = 0.
222  wake_s(:)        = 0.
223  wake_cstar(:)    = 0.
224  wake_fip(:)      = 0.
225  wake_pe          = 0.
226  fm_therm         = 0.
227  entr_therm       = 0.
228  detr_therm       = 0.
229
230  CALL fonte_neige_init(run_off_lic_0)
231  CALL pbl_surface_init( fder, snsrf, qsolsrf, tsoil )
232  CALL phyredem( "startphy.nc" )
233
234!  WRITE(lunout,*)'CCCCCCCCCCCCCCCCCC REACTIVER SORTIE VISU DANS ETAT0'
235!  WRITE(lunout,*)'entree histclo'
236  CALL histclo()
237
238END SUBROUTINE etat0phys_netcdf
239!
240!-------------------------------------------------------------------------------
241
242
243!-------------------------------------------------------------------------------
244!
245SUBROUTINE start_init_orog(lon_in,lat_in,phis,masque)
246!
247!===============================================================================
248! Comment:
249!   This routine launch grid_noro, which computes parameters for SSO scheme as
250!   described in LOTT & MILLER (1997) and LOTT(1999).
251!===============================================================================
252  USE grid_noro_m, ONLY: grid_noro
253  IMPLICIT NONE
254!-------------------------------------------------------------------------------
255! Arguments:
256  REAL,    INTENT(IN)    :: lon_in(:), lat_in(:)   ! dim (iml) (jml)
257  REAL,    INTENT(INOUT) :: phis(:,:), masque(:,:) ! dim (iml,jml)
258!-------------------------------------------------------------------------------
259! Local variables:
260  CHARACTER(LEN=256) :: modname
261  INTEGER            :: fid, llm_tmp,ttm_tmp, iml,jml, iml_rel,jml_rel, itau(1)
262  REAL               :: lev(1), date, dt
263  REAL, ALLOCATABLE  :: lon_rad(:), lon_ini(:), lon_rel(:,:), relief_hi(:,:)
264  REAL, ALLOCATABLE  :: lat_rad(:), lat_ini(:), lat_rel(:,:), tmp_var  (:,:)
265  REAL, ALLOCATABLE  :: zmea0(:,:), zstd0(:,:), zsig0(:,:)
266  REAL, ALLOCATABLE  :: zgam0(:,:), zthe0(:,:), zpic0(:,:), zval0(:,:)
267!-------------------------------------------------------------------------------
268  modname="start_init_orog"
269  iml=assert_eq(SIZE(lon_in),SIZE(phis,1),SIZE(masque,1),TRIM(modname)//" iml")
270  jml=assert_eq(SIZE(lat_in),SIZE(phis,2),SIZE(masque,2),TRIM(modname)//" jml")
271
272!--- HIGH RESOLUTION OROGRAPHY
273  CALL flininfo(orofname, iml_rel, jml_rel, llm_tmp, ttm_tmp, fid)
274
275  ALLOCATE(lat_rel(iml_rel,jml_rel),lon_rel(iml_rel,jml_rel))
276  CALL flinopen(orofname, .FALSE., iml_rel, jml_rel, llm_tmp, lon_rel, lat_rel,&
277                lev, ttm_tmp, itau, date, dt, fid)
278  ALLOCATE(relief_hi(iml_rel,jml_rel))
279  CALL flinget(fid, orogvar, iml_rel, jml_rel, llm_tmp, ttm_tmp, 1,1, relief_hi)
280  CALL flinclo(fid)
281
282!--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS
283  ALLOCATE(lon_ini(iml_rel),lat_ini(jml_rel))
284  lon_ini(:)=lon_rel(:,1); IF(MAXVAL(lon_rel)>pi) lon_ini=lon_ini*deg2rad
285  lat_ini(:)=lat_rel(1,:); IF(MAXVAL(lat_rel)>pi) lat_ini=lat_ini*deg2rad
286
287!--- FIELDS ARE PROCESSED TO BE ON STANDARD ANGULAR DOMAINS
288  ALLOCATE(lon_rad(iml_rel),lat_rad(jml_rel))
289  CALL conf_dat2d(orogvar, lon_ini, lat_ini, lon_rad, lat_rad, relief_hi,.FALSE.)
290  DEALLOCATE(lon_ini,lat_ini)
291
292!--- COMPUTING THE REQUIRED FIELDS USING ROUTINE grid_noro
293  WRITE(lunout,*)
294  WRITE(lunout,*)'*** Compute parameters needed for gravity wave drag code ***'
295
296!--- ALLOCATIONS OF SUB-CELL SCALES QUANTITIES
297  ALLOCATE(zmea0(iml,jml),zstd0(iml,jml)) !--- Mean orography and std deviation
298  ALLOCATE(zsig0(iml,jml),zgam0(iml,jml)) !--- Slope and nisotropy
299  ALLOCATE(zthe0(iml,jml))                !--- Highest slope orientation
300  ALLOCATE(zpic0(iml,jml),zval0(iml,jml)) !--- Peaks and valley heights
301
302!--- CALL OROGRAPHY MODULE TO COMPUTE FIELDS
303  CALL grid_noro(lon_rad,lat_rad,relief_hi,lon_in,lat_in,phis,zmea0,zstd0,     &
304                                      zsig0,zgam0,zthe0,zpic0,zval0,masque)
305  phis = phis * 9.81
306  phis(iml,:) = phis(1,:)
307  DEALLOCATE(relief_hi,lon_rad,lat_rad)
308
309!--- PUT QUANTITIES TO PHYSICAL GRID
310  CALL gr_dyn_fi(1,iml,jml,klon,zmea0,zmea); DEALLOCATE(zmea0)
311  CALL gr_dyn_fi(1,iml,jml,klon,zstd0,zstd); DEALLOCATE(zstd0)
312  CALL gr_dyn_fi(1,iml,jml,klon,zsig0,zsig); DEALLOCATE(zsig0)
313  CALL gr_dyn_fi(1,iml,jml,klon,zgam0,zgam); DEALLOCATE(zgam0)
314  CALL gr_dyn_fi(1,iml,jml,klon,zthe0,zthe); DEALLOCATE(zthe0)
315  CALL gr_dyn_fi(1,iml,jml,klon,zpic0,zpic); DEALLOCATE(zpic0)
316  CALL gr_dyn_fi(1,iml,jml,klon,zval0,zval); DEALLOCATE(zval0)
317
318
319END SUBROUTINE start_init_orog
320!
321!-------------------------------------------------------------------------------
322
323
324!-------------------------------------------------------------------------------
325!
326SUBROUTINE start_init_phys(lon_in,lat_in,phis)
327!
328!===============================================================================
329! Purpose:   Compute tsol and qsol, knowing phis.
330!===============================================================================
331  IMPLICIT NONE
332!-------------------------------------------------------------------------------
333! Arguments:
334  REAL,    INTENT(IN) :: lon_in(:),  lat_in(:)       ! dim (iml) (jml2)
335  REAL,    INTENT(IN) :: phis(:,:)                   ! dim (iml,jml)
336!-------------------------------------------------------------------------------
337! Local variables:
338  CHARACTER(LEN=256) :: modname
339  REAL               :: date, dt
340  INTEGER            :: iml, jml, jml2, itau(1)
341  REAL, ALLOCATABLE  :: lon_rad(:), lon_ini(:), var_ana(:,:)
342  REAL, ALLOCATABLE  :: lat_rad(:), lat_ini(:)
343  REAL, ALLOCATABLE  :: ts(:,:), qs(:,:)
344!-------------------------------------------------------------------------------
345  modname="start_init_phys"
346  iml=assert_eq(SIZE(lon_in),SIZE(phis,1),TRIM(modname)//" iml")
347  jml=SIZE(phis,2); jml2=SIZE(lat_in)
348
349  WRITE(lunout,*)'Opening the surface analysis'
350  CALL flininfo(phyfname, iml_phys, jml_phys, llm_phys, ttm_phys, fid_phys)
351  WRITE(lunout,*) 'Values read: ',  iml_phys, jml_phys, llm_phys, ttm_phys
352
353  ALLOCATE(lat_phys(iml_phys,jml_phys),lon_phys(iml_phys,jml_phys))
354  ALLOCATE(levphys_ini(llm_phys))
355  CALL flinopen(phyfname, .FALSE., iml_phys, jml_phys, llm_phys,              &
356                lon_phys,lat_phys,levphys_ini,ttm_phys,itau,date,dt,fid_phys)
357
358!--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS
359  ALLOCATE(lon_ini(iml_phys),lat_ini(jml_phys))
360  lon_ini(:)=lon_phys(:,1); IF(MAXVAL(lon_phys)>pi) lon_ini=lon_ini*deg2rad
361  lat_ini(:)=lat_phys(1,:); IF(MAXVAL(lat_phys)>pi) lat_ini=lat_ini*deg2rad
362
363  ALLOCATE(var_ana(iml_phys,jml_phys),lon_rad(iml_phys),lat_rad(jml_phys))
364  CALL get_var_phys(tsrfvar,ts)                   !--- SURFACE TEMPERATURE
365  CALL get_var_phys(qsolvar,qs)                   !--- SOIL MOISTURE
366  CALL flinclo(fid_phys)
367  DEALLOCATE(var_ana,lon_rad,lat_rad,lon_ini,lat_ini)
368
369!--- TSOL AND QSOL ON PHYSICAL GRID
370  ALLOCATE(tsol(klon))
371  CALL gr_dyn_fi(1,iml,jml,klon,ts,tsol)
372  CALL gr_dyn_fi(1,iml,jml,klon,qs,qsol)
373  DEALLOCATE(ts,qs)
374
375CONTAINS
376
377!-------------------------------------------------------------------------------
378!
379SUBROUTINE get_var_phys(title,field)
380!
381!-------------------------------------------------------------------------------
382  IMPLICIT NONE
383!-------------------------------------------------------------------------------
384! Arguments:
385  CHARACTER(LEN=*),  INTENT(IN)    :: title
386  REAL, ALLOCATABLE, INTENT(INOUT) :: field(:,:)
387!-------------------------------------------------------------------------------
388! Local variables:
389  INTEGER :: tllm
390!-------------------------------------------------------------------------------
391  SELECT CASE(title)
392    CASE(psrfvar);         tllm=0
393    CASE(tsrfvar,qsolvar); tllm=llm_phys
394  END SELECT
395  IF(ALLOCATED(field)) RETURN
396  ALLOCATE(field(iml,jml)); field(:,:)=0.
397  CALL flinget(fid_phys,title,iml_phys,jml_phys,tllm,ttm_phys,1,1,var_ana)
398  CALL conf_dat2d(title, lon_ini, lat_ini, lon_rad, lat_rad, var_ana, .TRUE.)
399  CALL interp_startvar(title, .TRUE., lon_rad, lat_rad, var_ana,               &
400                                      lon_in,  lat_in,  field)
401
402END SUBROUTINE get_var_phys
403!
404!-------------------------------------------------------------------------------
405!
406END SUBROUTINE start_init_phys
407!
408!-------------------------------------------------------------------------------
409
410
411!-------------------------------------------------------------------------------
412!
413SUBROUTINE interp_startvar(nam,ibeg,lon,lat,vari,lon2,lat2,varo)
414!
415!-------------------------------------------------------------------------------
416  USE inter_barxy_m, ONLY: inter_barxy
417  IMPLICIT NONE
418!-------------------------------------------------------------------------------
419! Arguments:
420  CHARACTER(LEN=*), INTENT(IN)  :: nam
421  LOGICAL,          INTENT(IN)  :: ibeg
422  REAL,             INTENT(IN)  :: lon(:), lat(:)   ! dim (ii) (jj)
423  REAL,             INTENT(IN)  :: vari(:,:)        ! dim (ii,jj)
424  REAL,             INTENT(IN)  :: lon2(:), lat2(:) ! dim (i1) (j2)
425  REAL,             INTENT(OUT) :: varo(:,:)        ! dim (i1) (j1)
426!-------------------------------------------------------------------------------
427! Local variables:
428  CHARACTER(LEN=256) :: modname
429  INTEGER            :: ii, jj, i1, j1, j2
430  REAL, ALLOCATABLE  :: vtmp(:,:)
431!-------------------------------------------------------------------------------
432  modname="interp_startvar"
433  ii=assert_eq(SIZE(lon), SIZE(vari,1),TRIM(modname)//" ii")
434  jj=assert_eq(SIZE(lat), SIZE(vari,2),TRIM(modname)//" jj")
435  i1=assert_eq(SIZE(lon2),SIZE(varo,1),TRIM(modname)//" i1")
436  j1=SIZE(varo,2); j2=SIZE(lat2)
437  ALLOCATE(vtmp(i1-1,j1))
438  IF(ibeg.AND.prt_level>1) THEN
439    WRITE(lunout,*)"--------------------------------------------------------"
440    WRITE(lunout,*)"$$$ Interpolation barycentrique pour "//TRIM(nam)//" $$$"
441    WRITE(lunout,*)"--------------------------------------------------------"
442  END IF
443  CALL inter_barxy(lon, lat(:jj-1), vari, lon2(:i1-1), lat2, vtmp)
444  CALL gr_int_dyn(vtmp, varo, i1-1, j1)
445
446END SUBROUTINE interp_startvar
447!
448!-------------------------------------------------------------------------------
449
450
451END MODULE etat0phys
452!
453!*******************************************************************************
454
Note: See TracBrowser for help on using the repository browser.