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

Last change on this file since 2320 was 2320, checked in by Ehouarn Millour, 9 years ago

Physics/dynamics separation: make an infotrac_phy module, which should be used from within the physics, and is initialized from infotrac (dynamics) via iniphysiq.
EM

File size: 22.6 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
39  USE phys_state_var_mod, ONLY: zmea, zstd, zsig, zgam, zthe, zpic, zval, z0m, &
40    rlon, solsw, radsol, t_ancien, wake_deltat, wake_s,  rain_fall, qsol, z0h, &
41    rlat, sollw, rugoro, q_ancien, wake_deltaq, wake_pe, snow_fall, ratqs,w01, &
42    sig1, ftsol, clwcon, fm_therm, wake_Cstar,  pctsrf,  entr_therm,radpas, f0,&
43    zmax0,fevap, rnebcon,falb_dir, wake_fip,    agesno,  detr_therm, pbl_tke,  &
44    phys_state_var_init
45
46  PRIVATE
47  PUBLIC :: etat0phys_netcdf
48
49  include "iniprint.h"
50  include "dimensions.h"
51  include "paramet.h"
52  include "comgeom2.h"
53  include "comvert.h"
54  include "comconst.h"
55  include "dimsoil.h"
56  include "temps.h"
57  include "comdissnew.h"
58  include "serre.h"
59  include "clesphys.h"
60  REAL, SAVE :: deg2rad
61  REAL, SAVE, ALLOCATABLE :: tsol(:)
62!  REAL, SAVE, ALLOCATABLE :: rugo(:,:)  ! ??? COMPUTED BUT NOT USED ???
63  INTEGER,            SAVE      :: iml_phys, jml_phys, llm_phys, ttm_phys, fid_phys
64  REAL, ALLOCATABLE,  SAVE      :: lon_phys(:,:), lat_phys(:,:), levphys_ini(:)
65  CHARACTER(LEN=256), PARAMETER :: orofname="Relief.nc"
66  CHARACTER(LEN=256), PARAMETER :: title="RELIEF"
67
68
69CONTAINS
70
71
72!-------------------------------------------------------------------------------
73!
74SUBROUTINE etat0phys_netcdf(ib, masque, phis)
75!
76!-------------------------------------------------------------------------------
77! Purpose: Creates initial states
78!-------------------------------------------------------------------------------
79! Note: This routine is designed to work for Earth
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 exner_hyb_m,    ONLY: exner_hyb
88  USE exner_milieu_m, ONLY: exner_milieu
89  USE test_disvert_m, ONLY: test_disvert
90  USE grid_atob_m,    ONLY: grille_m
91  IMPLICIT NONE
92!-------------------------------------------------------------------------------
93! Arguments:
94  LOGICAL, INTENT(IN)    :: ib          !--- Barycentric interpolation
95  REAL,    INTENT(INOUT) :: masque(:,:) !--- Land mask           dim(iip1,jjp1)
96  REAL,    INTENT(INOUT) :: phis  (:,:) !--- Ground geopotential dim(iip1,jjp1)
97!-------------------------------------------------------------------------------
98! Local variables:
99  CHARACTER(LEN=256) :: modname="etat0phys_netcdf", fmt
100  INTEGER            :: i, j, l, ji, iml, jml
101  REAL               :: phystep
102  REAL, DIMENSION(SIZE(masque,1),SIZE(masque,2)) :: masque_tmp
103  REAL, DIMENSION(klon)               :: sn, rugmer, run_off_lic_0, fder
104  REAL, DIMENSION(klon,nbsrf)         :: qsolsrf, snsrf
105  REAL, DIMENSION(klon,nsoilmx,nbsrf) :: tsoil
106
107!--- Local variables for sea-ice reading:
108  LOGICAL           :: read_mask
109  INTEGER           :: iml_lic, jml_lic, isst(klon-2)
110  INTEGER           :: fid, llm_tmp, ttm_tmp, itaul(1)
111  REAL, ALLOCATABLE :: dlon_lic(:), lon_lic(:,:), fraclic(:,:)
112  REAL, ALLOCATABLE :: dlat_lic(:), lat_lic(:,:)
113  REAL              :: date, lev(1), dummy
114  REAL              :: flic_tmp(SIZE(masque,1),SIZE(masque,2))
115
116!--- Arguments for conf_phys
117  LOGICAL :: ok_journe, ok_mensuel, ok_instan, ok_hf, ok_LES, callstats
118  REAL    :: solarlong0, seuil_inversion, fact_cldcon, facttemps
119  LOGICAL :: ok_newmicro
120  INTEGER :: iflag_radia, iflag_cldcon, iflag_ratqs
121  REAL    :: ratqsbas, ratqshaut, tau_ratqs
122  LOGICAL :: ok_ade, ok_aie, ok_cdnc, aerosol_couple
123  INTEGER :: flag_aerosol
124  LOGICAL :: flag_aerosol_strat
125  LOGICAL :: new_aod
126  REAL    :: bl95_b0, bl95_b1
127  INTEGER :: read_climoz                        !--- Read ozone climatology
128  REAL    :: alp_offset
129
130  deg2rad= pi/180.0
131  iml=assert_eq(SIZE(masque,1),SIZE(phis,1),TRIM(modname)//" iml")
132  jml=assert_eq(SIZE(masque,2),SIZE(phis,2),TRIM(modname)//" jml")
133
134! Grid construction and miscellanous initializations.
135!*******************************************************************************
136  CALL conf_phys(  ok_journe, ok_mensuel, ok_instan, ok_hf, ok_LES,     &
137                   callstats,                                           &
138                   solarlong0,seuil_inversion,                          &
139                   fact_cldcon, facttemps,ok_newmicro,iflag_radia,      &
140                   iflag_cldcon,                                        &
141                   iflag_ratqs,ratqsbas,ratqshaut,tau_ratqs,            &
142                   ok_ade, ok_aie, ok_cdnc, aerosol_couple,             &
143                   flag_aerosol, flag_aerosol_strat, new_aod,           &
144                   bl95_b0, bl95_b1,                                    &
145                   read_climoz,                                         &
146                   alp_offset)
147
148  CALL phys_state_var_init(read_climoz)
149
150  co2_ppm0 = co2_ppm  !--- Initial atmospheric CO2 conc. from .def file
151
152  rlat(1) = pi/2.
153  DO j=2,jjm; rlat((j-2)*iim+2:(j-1)*iim+1)=rlatu(j);     END DO
154  rlat(klon) = - pi/2.
155  rlat(:)=rlat(:)*(180.0/pi)
156
157  rlon(1) = 0.0
158  DO j=2,jjm; rlon((j-2)*iim+2:(j-1)*iim+1)=rlonv(1:iim); END DO
159  rlon(klon) = 0.0
160  rlon(:)=rlon(:)*(180.0/pi)
161
162! Compute ground geopotential, sub-cells quantities and possibly the mask.
163!*******************************************************************************
164  read_mask=ANY(masque/=-99999.); masque_tmp=masque
165  CALL start_init_orog(rlonv, rlatu, phis, masque_tmp)
166  WRITE(fmt,"(i4,'i1)')")iml ; fmt='('//ADJUSTL(fmt)
167  IF(.NOT.read_mask) THEN                       !--- Keep mask form orography
168    masque=masque_tmp
169    IF(prt_level>=1) THEN
170      WRITE(lunout,*)'BUILT MASK :'
171      WRITE(lunout,fmt) NINT(masque)
172    END IF
173    WHERE(   masque(:,:)<EPSFRA) masque(:,:)=0.
174    WHERE(1.-masque(:,:)<EPSFRA) masque(:,:)=1.
175  END IF
176 CALL gr_dyn_fi(1,iml,jml,klon,masque,zmasq) !--- Land mask to physical grid
177
178! Compute tsol and qsol on physical grid, knowing phis on 2D grid.
179!*******************************************************************************
180  CALL start_init_phys(rlonv, rlatu, rlonu, rlatv, ib, phis)
181
182! Some initializations.
183!*******************************************************************************
184  sn    (:) = 0.0                               !--- Snow
185  radsol(:) = 0.0                               !--- Net radiation at ground
186  rugmer(:) = 0.001                             !--- Ocean rugosity
187  IF(read_climoz>=1) &                          !--- Ozone climatology
188    CALL regr_lat_time_climoz(read_climoz)
189
190! Sub-surfaces initialization
191!*******************************************************************************
192!--- Read and interpolate on model T-grid soil fraction and soil ice fraction.
193  CALL flininfo("landiceref.nc", iml_lic, jml_lic, llm_tmp, ttm_tmp, fid)
194  ALLOCATE( lat_lic(iml_lic,jml_lic),lon_lic(iml_lic, jml_lic))
195  ALLOCATE(dlat_lic(jml_lic),       dlon_lic(iml_lic))
196  ALLOCATE( fraclic(iml_lic,jml_lic))
197  CALL flinopen("landiceref.nc", .FALSE., iml_lic, jml_lic, llm_tmp,  &
198 &               lon_lic, lat_lic, lev, ttm_tmp, itaul, date, dt, fid)
199  CALL flinget(fid, 'landice', iml_lic, jml_lic, llm_tmp, ttm_tmp, 1,1, fraclic)
200  CALL flinclo(fid)
201  WRITE(lunout,*)'landice dimensions: iml_lic, jml_lic : ',iml_lic,jml_lic
202  IF(MAXVAL(lon_lic)>pi) lon_lic=lon_lic*pi/180. !--- Conversion to degrees
203  IF(MAXVAL(lat_lic)>pi) lat_lic=lat_lic*pi/180.
204  dlon_lic(:)=lon_lic(:,1)
205  dlat_lic(:)=lat_lic(1,:)
206  CALL grille_m(dlon_lic, dlat_lic, fraclic, rlonv(1:iim), rlatu, flic_tmp(1:iim,:) )
207  flic_tmp(iml,:)=flic_tmp(1,:)
208
209!--- To the physical grid
210  pctsrf(:,:) = 0.
211  CALL gr_dyn_fi(1, iml, jml, klon, flic_tmp, pctsrf(:,is_lic))
212
213!--- Adequation with soil/sea mask
214  WHERE(pctsrf(:,is_lic)<EPSFRA) pctsrf(:,is_lic)=0.
215  WHERE(zmasq(:)<EPSFRA)         pctsrf(:,is_lic)=0.
216  pctsrf(:,is_ter)=zmasq(:)
217  DO ji=1,klon
218    IF(zmasq(ji)>EPSFRA) THEN
219      IF(pctsrf(ji,is_lic)>=zmasq(ji)) THEN
220        pctsrf(ji,is_lic)=zmasq(ji)
221        pctsrf(ji,is_ter)=0.
222      ELSE
223        pctsrf(ji,is_ter)=zmasq(ji)-pctsrf(ji,is_lic)
224        IF(pctsrf(ji,is_ter)<EPSFRA) THEN
225          pctsrf(ji,is_ter)=0.
226          pctsrf(ji,is_lic)=zmasq(ji)
227        END IF
228      END IF
229    END IF
230  END DO
231
232!--- Sub-surface ocean and sea ice (sea ice set to zero for start).
233  pctsrf(:,is_oce)=(1.-zmasq(:))
234  WHERE(pctsrf(:,is_oce)<EPSFRA) pctsrf(:,is_oce)=0.
235  IF(read_mask) pctsrf(:,is_oce)=1-zmasq(:)
236  isst=0
237  WHERE(pctsrf(2:klon-1,is_oce)>0.) isst=1
238
239!--- It is checked that the sub-surfaces sum is equal to 1.
240  ji=COUNT((ABS(SUM(pctsrf(:,:),dim=2))-1.0)>EPSFRA)
241  IF(ji/=0) WRITE(lunout,*) 'Sub-cell distribution problem for ',ji,' points'
242
243! Write physical initial state
244!*******************************************************************************
245  WRITE(lunout,*)'phystep ',dtvr,iphysiq,nbapp_rad
246  phystep = dtvr * FLOAT(iphysiq)
247  radpas  = NINT (86400./phystep/ FLOAT(nbapp_rad) )
248  WRITE(lunout,*)'phystep =', phystep, radpas
249
250! Init: ftsol, snsrf, qsolsrf, tsoil, rain_fall, snow_fall, solsw, sollw, z0
251!*******************************************************************************
252  DO i=1,nbsrf; ftsol(:,i) = tsol; END DO
253  DO i=1,nbsrf; snsrf(:,i) = sn;   END DO
254  falb_dir(:,is_ter,:) = 0.08
255  falb_dir(:,is_lic,:) = 0.6
256  falb_dir(:,is_oce,:) = 0.5
257  falb_dir(:,is_sic,:) = 0.6
258  fevap(:,:) = 0.
259  DO i=1,nbsrf; qsolsrf(:,i)=150.; END DO
260  DO i=1,nbsrf; DO j=1,nsoilmx; tsoil(:,j,i) = tsol; END DO; END DO
261  rain_fall  = 0.
262  snow_fall  = 0.
263  solsw      = 165.
264  sollw      = -53.
265  t_ancien   = 273.15
266  q_ancien   = 0.
267  agesno     = 0.
268
269  z0m(:,is_oce) = rugmer(:)
270  z0m(:,is_ter) = MAX(1.0e-05,zstd(:)*zsig(:)/2.0)
271  z0m(:,is_lic) = MAX(1.0e-05,zstd(:)*zsig(:)/2.0)
272  z0m(:,is_sic) = 0.001
273  z0h(:,:)=z0m(:,:)
274
275  fder    = 0.0
276  clwcon  = 0.0
277  rnebcon = 0.0
278  ratqs   = 0.0
279  run_off_lic_0 = 0.0
280  rugoro  = 0.0
281
282! Before phyredem calling, surface modules and values to be saved in startphy.nc
283! are initialized
284!*******************************************************************************
285  dummy            = 1.0
286  pbl_tke(:,:,:)   = 1.e-8
287  zmax0(:)         = 40.
288  f0(:)            = 1.e-5
289  sig1(:,:)        = 0.
290  w01(:,:)         = 0.
291  wake_deltat(:,:) = 0.
292  wake_deltaq(:,:) = 0.
293  wake_s(:)        = 0.
294  wake_cstar(:)    = 0.
295  wake_fip(:)      = 0.
296  wake_pe          = 0.
297  fm_therm         = 0.
298  entr_therm       = 0.
299  detr_therm       = 0.
300
301  CALL fonte_neige_init(run_off_lic_0)
302  CALL pbl_surface_init( fder, snsrf, qsolsrf, tsoil )
303  CALL phyredem( "startphy.nc" )
304
305!  WRITE(lunout,*)'CCCCCCCCCCCCCCCCCC REACTIVER SORTIE VISU DANS ETAT0'
306!  WRITE(lunout,*)'entree histclo'
307  CALL histclo()
308
309END SUBROUTINE etat0phys_netcdf
310!
311!-------------------------------------------------------------------------------
312
313
314!-------------------------------------------------------------------------------
315!
316SUBROUTINE start_init_orog(lon_in,lat_in,phis,masque)
317!
318!===============================================================================
319! Comment:
320!   This routine launch grid_noro, which computes parameters for SSO scheme as
321!   described in LOTT & MILLER (1997) and LOTT(1999).
322!===============================================================================
323  USE conf_dat_m,  ONLY: conf_dat2d
324!  USE grid_atob_m, ONLY: rugsoro
325  USE grid_noro_m, ONLY: grid_noro
326  IMPLICIT NONE
327!-------------------------------------------------------------------------------
328! Arguments:
329  REAL,    INTENT(IN)    :: lon_in(:), lat_in(:)   ! dim (iml) (jml)
330  REAL,    INTENT(INOUT) :: phis(:,:), masque(:,:) ! dim (iml,jml)
331!-------------------------------------------------------------------------------
332! Local variables:
333  CHARACTER(LEN=256) :: modname="start_init_orog"
334  CHARACTER(LEN=256) :: title="RELIEF"
335  INTEGER            :: fid, llm_tmp,ttm_tmp, iml,jml, iml_rel,jml_rel, itau(1)
336  REAL               :: lev(1), date, dt
337  REAL, ALLOCATABLE  :: lon_rad(:), lon_ini(:), lon_rel(:,:), relief_hi(:,:)
338  REAL, ALLOCATABLE  :: lat_rad(:), lat_ini(:), lat_rel(:,:), tmp_var  (:,:)
339  REAL, ALLOCATABLE  :: zmea0(:,:), zstd0(:,:), zsig0(:,:)
340  REAL, ALLOCATABLE  :: zgam0(:,:), zthe0(:,:), zpic0(:,:), zval0(:,:)
341!-------------------------------------------------------------------------------
342  iml=assert_eq(SIZE(lon_in),SIZE(phis,1),SIZE(masque,1),TRIM(modname)//" iml")
343  jml=assert_eq(SIZE(lat_in),SIZE(phis,2),SIZE(masque,2),TRIM(modname)//" jml")
344
345!--- HIGH RESOLUTION OROGRAPHY
346  CALL flininfo(orofname, iml_rel, jml_rel, llm_tmp, ttm_tmp, fid)
347
348  ALLOCATE(lat_rel(iml_rel,jml_rel),lon_rel(iml_rel,jml_rel))
349  CALL flinopen(orofname, .FALSE., iml_rel, jml_rel, llm_tmp, lon_rel, lat_rel,&
350                lev, ttm_tmp, itau, date, dt, fid)
351  ALLOCATE(relief_hi(iml_rel,jml_rel))
352  CALL flinget(fid, title, iml_rel, jml_rel, llm_tmp, ttm_tmp, 1, 1, relief_hi)
353  CALL flinclo(fid)
354
355!--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS
356  ALLOCATE(lon_ini(iml_rel),lat_ini(jml_rel))
357  lon_ini(:)=lon_rel(:,1); IF(MAXVAL(lon_rel)>pi) lon_ini=lon_ini*deg2rad
358  lat_ini(:)=lat_rel(1,:); IF(MAXVAL(lat_rel)>pi) lat_ini=lat_ini*deg2rad
359
360!--- FIELDS ARE PROCESSED TO BE ON STANDARD ANGULAR DOMAINS
361  ALLOCATE(lon_rad(iml_rel),lat_rad(jml_rel))
362  CALL conf_dat2d(title, lon_ini, lat_ini, lon_rad, lat_rad, relief_hi, .FALSE.)
363  DEALLOCATE(lon_ini,lat_ini)
364
365!--- COMPUTING THE REQUIRED FIELDS USING ROUTINE grid_noro
366  WRITE(lunout,*)
367  WRITE(lunout,*)'*** Compute parameters needed for gravity wave drag code ***'
368
369!--- ALLOCATIONS OF SUB-CELL SCALES QUANTITIES
370  ALLOCATE(zmea0(iml,jml),zstd0(iml,jml)) !--- Mean orography and std deviation
371  ALLOCATE(zsig0(iml,jml),zgam0(iml,jml)) !--- Slope and nisotropy
372  ALLOCATE(zthe0(iml,jml))                !--- Highest slope orientation
373  ALLOCATE(zpic0(iml,jml),zval0(iml,jml)) !--- Peaks and valley heights
374
375!--- CALL OROGRAPHY MODULE TO COMPUTE FIELDS
376  CALL grid_noro(lon_rad,lat_rad,relief_hi,lon_in,lat_in,phis,zmea0,zstd0,     &
377                                      zsig0,zgam0,zthe0,zpic0,zval0,masque)
378  phis = phis * 9.81
379  phis(iml,:) = phis(1,:)
380
381!--- COMPUTE SURFACE ROUGHNESS
382!  WRITE(lunout,*)
383!  WRITE(lunout,*)'*** Compute surface roughness induced by the orography ***'
384!  ALLOCATE(tmp_var(iml-1,jml))
385!  CALL rugsoro(lon_rad, lat_rad, relief_hi, lon_in(1:iml-1), lat_in, tmp_var)
386!  ALLOCATE(rugo(iml,jml)); rugo(1:iml-1,:)=tmp_var; rugo(iml,:)=tmp_var(1,:)
387!  DEALLOCATE(tmp_var)
388  DEALLOCATE(relief_hi,lon_rad,lat_rad)
389
390!--- PUT QUANTITIES TO PHYSICAL GRID
391  CALL gr_dyn_fi(1,iml,jml,klon,zmea0,zmea); DEALLOCATE(zmea0)
392  CALL gr_dyn_fi(1,iml,jml,klon,zstd0,zstd); DEALLOCATE(zstd0)
393  CALL gr_dyn_fi(1,iml,jml,klon,zsig0,zsig); DEALLOCATE(zsig0)
394  CALL gr_dyn_fi(1,iml,jml,klon,zgam0,zgam); DEALLOCATE(zgam0)
395  CALL gr_dyn_fi(1,iml,jml,klon,zthe0,zthe); DEALLOCATE(zthe0)
396  CALL gr_dyn_fi(1,iml,jml,klon,zpic0,zpic); DEALLOCATE(zpic0)
397  CALL gr_dyn_fi(1,iml,jml,klon,zval0,zval); DEALLOCATE(zval0)
398
399
400END SUBROUTINE start_init_orog
401!
402!-------------------------------------------------------------------------------
403
404
405!-------------------------------------------------------------------------------
406!
407SUBROUTINE start_init_phys(lon_in,lat_in,lon_in2,lat_in2,ibar,phis)
408!
409!===============================================================================
410! Purpose:   Compute tsol and qsol, knowing phis.
411!===============================================================================
412  IMPLICIT NONE
413!-------------------------------------------------------------------------------
414! Arguments:
415  REAL,    INTENT(IN) :: lon_in (:),  lat_in (:)     ! dim (iml) (jml)
416  REAL,    INTENT(IN) :: lon_in2(:),  lat_in2(:)     ! dim (iml) (jml2)
417  LOGICAL, INTENT(IN) :: ibar
418  REAL,    INTENT(IN) :: phis(:,:)                   ! dim (iml,jml)
419!-------------------------------------------------------------------------------
420! Local variables:
421  CHARACTER(LEN=256) :: modname="start_init_phys", physfname="ECPHY.nc"
422  REAL               :: date, dt
423  INTEGER            :: iml, jml, jml2, itau(1)
424  REAL, ALLOCATABLE  :: lon_rad(:), lon_ini(:), var_ana(:,:)
425  REAL, ALLOCATABLE  :: lat_rad(:), lat_ini(:)
426  REAL, ALLOCATABLE  :: ts(:,:), qs(:,:)
427!-------------------------------------------------------------------------------
428  iml=assert_eq(SIZE(lon_in),SIZE(phis,1),SIZE(lon_in2),TRIM(modname)//" iml")
429  jml=assert_eq(SIZE(lat_in),SIZE(phis,2),              TRIM(modname)//" jml")
430  jml2=SIZE(lat_in2)
431
432  WRITE(lunout,*)'Opening the surface analysis'
433  CALL flininfo(physfname, iml_phys, jml_phys, llm_phys, ttm_phys, fid_phys)
434  WRITE(lunout,*) 'Values read: ',   iml_phys, jml_phys, llm_phys, ttm_phys
435
436  ALLOCATE(lat_phys(iml_phys,jml_phys),lon_phys(iml_phys,jml_phys))
437  ALLOCATE(levphys_ini(llm_phys))
438  CALL flinopen(physfname, .FALSE., iml_phys, jml_phys, llm_phys,              &
439                lon_phys,lat_phys,levphys_ini,ttm_phys,itau,date,dt,fid_phys)
440
441!--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS
442  ALLOCATE(lon_ini(iml_phys),lat_ini(jml_phys))
443  lon_ini(:)=lon_phys(:,1); IF(MAXVAL(lon_phys)>pi) lon_ini=lon_ini*deg2rad
444  lat_ini(:)=lat_phys(1,:); IF(MAXVAL(lat_phys)>pi) lat_ini=lat_ini*deg2rad
445
446  ALLOCATE(var_ana(iml_phys,jml_phys),lon_rad(iml_phys),lat_rad(jml_phys))
447  CALL get_var_phys('ST'  ,ts)                   !--- SURFACE TEMPERATURE
448  CALL get_var_phys('CDSW',qs)                   !--- SOIL MOISTURE
449  CALL flinclo(fid_phys)
450  DEALLOCATE(var_ana,lon_rad,lat_rad,lon_ini,lat_ini)
451
452!--- TSOL AND QSOL ON PHYSICAL GRID
453  ALLOCATE(tsol(klon))
454  CALL gr_dyn_fi(1,iml,jml,klon,ts,tsol)
455  CALL gr_dyn_fi(1,iml,jml,klon,qs,qsol)
456  DEALLOCATE(ts,qs)
457
458CONTAINS
459
460!-------------------------------------------------------------------------------
461!
462SUBROUTINE get_var_phys(title,field)
463!
464!-------------------------------------------------------------------------------
465  USE conf_dat_m, ONLY: conf_dat2d
466  IMPLICIT NONE
467!-------------------------------------------------------------------------------
468! Arguments:
469  CHARACTER(LEN=*),  INTENT(IN)    :: title
470  REAL, ALLOCATABLE, INTENT(INOUT) :: field(:,:)
471!-------------------------------------------------------------------------------
472! Local variables:
473  INTEGER :: tllm
474!-------------------------------------------------------------------------------
475  SELECT CASE(title)
476    CASE('SP');        tllm=0
477    CASE('ST','CDSW'); tllm=llm_phys
478  END SELECT
479  IF(ALLOCATED(field)) RETURN
480  ALLOCATE(field(iml,jml)); field(:,:)=0.
481  CALL flinget(fid_phys,title,iml_phys,jml_phys,tllm,ttm_phys,1,1,var_ana)
482  CALL conf_dat2d(title, lon_ini, lat_ini,  lon_rad, lat_rad, var_ana, ibar)
483  CALL interp_startvar(title, ibar, .TRUE., lon_rad, lat_rad, var_ana,         &
484                            lon_in, lat_in, lon_in2, lat_in2, field)
485
486END SUBROUTINE get_var_phys
487!
488!-------------------------------------------------------------------------------
489!
490END SUBROUTINE start_init_phys
491!
492!-------------------------------------------------------------------------------
493
494
495!-------------------------------------------------------------------------------
496!
497SUBROUTINE interp_startvar(nam,ibar,ibeg,lon,lat,vari,lon1,lat1,lon2,lat2,varo)
498!
499!-------------------------------------------------------------------------------
500  USE inter_barxy_m, ONLY: inter_barxy
501  USE grid_atob_m,   ONLY: grille_m
502  IMPLICIT NONE
503!-------------------------------------------------------------------------------
504! Arguments:
505  CHARACTER(LEN=*), INTENT(IN)  :: nam
506  LOGICAL,          INTENT(IN)  :: ibar, ibeg
507  REAL,             INTENT(IN)  :: lon(:), lat(:)   ! dim (ii) (jj)
508  REAL,             INTENT(IN)  :: vari(:,:)        ! dim (ii,jj)
509  REAL,             INTENT(IN)  :: lon1(:), lat1(:) ! dim (i1) (j1)
510  REAL,             INTENT(IN)  :: lon2(:), lat2(:) ! dim (i1) (j2)
511  REAL,             INTENT(OUT) :: varo(:,:)        ! dim (i1) (j1)
512!-------------------------------------------------------------------------------
513! Local variables:
514  CHARACTER(LEN=256) :: modname="interp_startvar"
515  INTEGER            :: ii, jj, i1, j1, j2
516  REAL, ALLOCATABLE  :: vtmp(:,:)
517!-------------------------------------------------------------------------------
518  ii=assert_eq(SIZE(lon),            SIZE(vari,1),TRIM(modname)//" ii")
519  jj=assert_eq(SIZE(lat),            SIZE(vari,2),TRIM(modname)//" jj")
520  i1=assert_eq(SIZE(lon1),SIZE(lon2),SIZE(varo,1),TRIM(modname)//" i1")
521  j1=assert_eq(SIZE(lat1),           SIZE(varo,2),TRIM(modname)//" j1")
522  j2=SIZE(lat2)
523  ALLOCATE(vtmp(i1-1,j1))
524  IF(ibar) THEN
525    IF(ibeg.AND.prt_level>1) THEN
526      WRITE(lunout,*)"--------------------------------------------------------"
527      WRITE(lunout,*)"$$$ Interpolation barycentrique pour "//TRIM(nam)//" $$$"
528      WRITE(lunout,*)"--------------------------------------------------------"
529    END IF
530    CALL inter_barxy(lon, lat(:jj-1), vari, lon2(:i1-1), lat2, vtmp)
531  ELSE
532    CALL grille_m   (lon, lat,        vari, lon1,        lat1, vtmp)
533  END IF
534  CALL gr_int_dyn(vtmp, varo, i1-1, j1)
535
536END SUBROUTINE interp_startvar
537!
538!-------------------------------------------------------------------------------
539
540
541END MODULE etat0phys
542!
543!*******************************************************************************
544
Note: See TracBrowser for help on using the repository browser.