source: LMDZ6/branches/Ocean_skin/libf/dynphy_lonlat/phylmd/etat0phys_netcdf.F90 @ 3605

Last change on this file since 3605 was 3605, checked in by lguez, 5 years ago

Merge revisions 3427:3600 of trunk into branch Ocean_skin

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