source: LMDZ6/trunk/libf/dynphy_lonlat/phylmd/etat0phys_netcdf.F90 @ 3520

Last change on this file since 3520 was 3520, checked in by oboucher, 5 years ago

Looks like ok_volcan was inserted in the wrong location

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