source: LMDZ5/trunk/libf/dynphy_lonlat/phylmd/etat0phys_netcdf.F90 @ 5382

Last change on this file since 5382 was 2899, checked in by lguez, 8 years ago

startphy.nc created by ce0l and lmdz1d was not reproducible because
ql_ancien, qs_ancien, prlw_ancien, prsw_ancien, prw_ancien were not
defined before the call to phyredem.

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