source: trunk/ICOSA_LMDZ/src/phystd/interface_icosa_lmdz.f90 @ 1700

Last change on this file since 1700 was 1698, checked in by emillour, 9 years ago

Dynamico-LMDZ interface:
Fix pseudoalt computation in interface_icosa_lmdz.f90
EM

File size: 27.0 KB
Line 
1MODULE interface_icosa_lmdz_mod
2
3  USE field_mod, ONLY: t_field
4  USE transfert_mod, ONLY: t_message
5 
6 
7  TYPE(t_message),SAVE :: req_u
8  TYPE(t_message),SAVE :: req_dps0, req_dulon0, req_dulat0, req_dTemp0, req_dq0
9
10  TYPE(t_field),POINTER,SAVE :: f_p(:)
11  TYPE(t_field),POINTER,SAVE :: f_pks(:) 
12  TYPE(t_field),POINTER,SAVE :: f_pk(:) 
13  TYPE(t_field),POINTER,SAVE :: f_p_layer(:)   
14  TYPE(t_field),POINTER,SAVE :: f_theta(:)   
15  TYPE(t_field),POINTER,SAVE :: f_phi(:)   
16  TYPE(t_field),POINTER,SAVE :: f_Temp(:)   
17  TYPE(t_field),POINTER,SAVE :: f_ulon(:)   
18  TYPE(t_field),POINTER,SAVE :: f_ulat(:)   
19  TYPE(t_field),POINTER,SAVE :: f_dulon(:)
20  TYPE(t_field),POINTER,SAVE :: f_dulat(:)
21  TYPE(t_field),POINTER,SAVE :: f_dTemp(:)
22  TYPE(t_field),POINTER,SAVE :: f_dq(:)
23  TYPE(t_field),POINTER,SAVE :: f_dps(:)
24  TYPE(t_field),POINTER,SAVE :: f_duc(:)
25  TYPE(t_field),POINTER,SAVE :: f_bounds_lon(:)
26  TYPE(t_field),POINTER,SAVE :: f_bounds_lat(:)
27
28  INTEGER,SAVE :: start_clock
29  INTEGER,SAVE :: stop_clock
30  INTEGER,SAVE :: count_clock=0
31 
32  REAL,SAVE :: day_length ! length of a day (s)
33  REAL,SAVE :: year_length ! length of a year (s)
34  INTEGER,SAVE :: start_day ! reference sol value at beginning of the run wrt Ls=0
35 
36  INTEGER,SAVE :: nbp_phys
37  INTEGER,SAVE :: nbp_phys_glo
38 
39  CHARACTER(len=20),SAVE,ALLOCATABLE :: tname(:) ! tracer names
40  REAL,SAVE :: pday ! number of ellapsed sols since Ls=0
41  REAL,SAVE :: ptime ! "universal time" as fraction of sol (e.g. 0.5 for noon)
42
43
44CONTAINS
45
46  SUBROUTINE initialize_physics
47  USE distrib_icosa_lmdz_mod, ONLY : init_distrib_icosa_lmdz, transfer_icosa_to_lmdz
48! from dynamico
49  USE domain_mod
50  USE dimensions
51  USE mpi_mod
52  USE mpipara
53  USE disvert_mod
54  USE xios_mod
55  USE time_mod , init_time_icosa=> init_time
56  USE transfert_mod
57 
58! from LMDZ
59  USE mod_grid_phy_lmdz, ONLY : unstructured
60  USE mod_phys_lmdz_omp_data, ONLY: klon_omp
61!  USE time_phylmdz_mod, ONLY: init_time_lmdz => init_time
62  USE transfert_mod
63  USE physics_distribution_mod, ONLY : init_physics_distribution
64!  USE geometry_mod, ONLY : init_geometry
65!  USE vertical_layers_mod, ONLY : init_vertical_layers
66!  USE infotrac_phy, ONLY : init_infotrac_phy
67!  USE inifis_mod, ONLY : inifis
68!  USE phyaqua_mod, ONLY : iniaqua
69   
70 
71  IMPLICIT NONE
72  INTEGER  :: ind,i,j,ij,pos
73  REAL(rstd),POINTER :: bounds_lon(:,:)
74  REAL(rstd),POINTER :: bounds_lat(:,:)
75 
76  REAL(rstd),ALLOCATABLE :: latfi(:)
77  REAL(rstd),ALLOCATABLE :: lonfi(:)
78  REAL(rstd),ALLOCATABLE :: airefi(:)
79  REAL(rstd),ALLOCATABLE :: bounds_latfi(:,:)
80  REAL(rstd),ALLOCATABLE :: bounds_lonfi(:,:)
81  REAL(rstd) :: pseudoalt(llm)
82
83  INTEGER :: run_length 
84  INTEGER :: annee_ref 
85  INTEGER :: day_ref   
86  INTEGER :: day_ini   
87  REAL    :: start_time
88  REAL    :: physics_timestep   
89
90
91  INTEGER                       :: nqo, nbtr
92  CHARACTER(len=4)              :: type_trac
93  CHARACTER(len=20),ALLOCATABLE :: tname(:)    ! tracer short name for restart and diagnostics
94  CHARACTER(len=23),ALLOCATABLE :: ttext(:)     ! tracer long name for diagnostics
95  INTEGER,ALLOCATABLE           :: niadv(:)    ! equivalent dyn / physique
96  INTEGER,ALLOCATABLE           :: conv_flg(:) ! conv_flg(it)=0 : convection desactivated for tracer number it
97  INTEGER,ALLOCATABLE           :: pbl_flg(:)  ! pbl_flg(it)=0  : boundary layer diffusion desactivaded for tracer number it
98  CHARACTER(len=8),ALLOCATABLE :: solsym(:)  ! tracer name from inca
99
100  INTEGER :: iflag_phys   
101  INTEGER :: nbp_phys, nbp_phys_glo
102 
103!$OMP PARALLEL
104    CALL allocate_field(f_bounds_lon,field_t,type_real,6)
105    CALL allocate_field(f_bounds_lat,field_t,type_real,6)
106    CALL allocate_field(f_p,field_t,type_real,llm+1,name="p_in")
107    CALL allocate_field(f_pks,field_t,type_real)
108    CALL allocate_field(f_pk,field_t,type_real,llm)
109    CALL allocate_field(f_p_layer,field_t,type_real,llm,name="p_layer_in")
110    CALL allocate_field(f_theta,field_t,type_real,llm)
111    CALL allocate_field(f_phi,field_t,type_real,llm,name="phi_in")
112    CALL allocate_field(f_Temp,field_t,type_real,llm,name="Temp_in")
113    CALL allocate_field(f_ulon,field_t,type_real,llm,name="ulon_in")
114    CALL allocate_field(f_ulat,field_t,type_real,llm,name="ulat_in")
115    CALL allocate_field(f_dulon,field_t,type_real,llm,name="dulon_out")
116    CALL allocate_field(f_dulat,field_t,type_real,llm,name="dulat_out")
117    CALL allocate_field(f_dTemp,field_t,type_real,llm,name="dTemp_out")
118    CALL allocate_field(f_dq,field_t,type_real,llm,nqtot,name="dq_out")
119    CALL allocate_field(f_dps,field_t,type_real,name="dps_out")
120    CALL allocate_field(f_duc,field_t,type_real,3,llm)   
121
122    CALL init_message(f_dps,req_i0,req_dps0)
123    CALL init_message(f_dulon,req_i0,req_dulon0)
124    CALL init_message(f_dulat,req_i0,req_dulat0)
125    CALL init_message(f_dTemp,req_i0,req_dTemp0)
126    CALL init_message(f_dq,req_i0,req_dq0)
127!$OMP END PARALLEL   
128
129    nbp_phys=0
130    DO ind=1,ndomain
131      CALL swap_dimensions(ind)
132      DO j=jj_begin,jj_end
133        DO i=ii_begin,ii_end
134          IF (domain(ind)%own(i,j)) nbp_phys=nbp_phys+1
135        ENDDO
136      ENDDO
137    ENDDO
138   
139
140!initialize LMDZ5 physic mpi decomposition
141    CALL MPI_ALLREDUCE(nbp_phys,nbp_phys_glo,1,MPI_INTEGER,MPI_SUM,comm_icosa,ierr)
142    CALL init_physics_distribution(unstructured, 6, nbp_phys, 1, nbp_phys_glo, llm, comm_icosa)
143   
144    DO ind=1,ndomain
145        CALL swap_dimensions(ind)
146        CALL swap_geometry(ind)
147        bounds_lon=f_bounds_lon(ind)
148        bounds_lat=f_bounds_lat(ind)
149        DO j=jj_begin,jj_end
150          DO i=ii_begin,ii_end
151            ij=(j-1)*iim+i
152            CALL xyz2lonlat(xyz_v(ij+z_rup,:), bounds_lon(ij,1), bounds_lat(ij,1))
153            CALL xyz2lonlat(xyz_v(ij+z_up,:), bounds_lon(ij,2), bounds_lat(ij,2))
154            CALL xyz2lonlat(xyz_v(ij+z_lup,:), bounds_lon(ij,3), bounds_lat(ij,3))
155            CALL xyz2lonlat(xyz_v(ij+z_ldown,:), bounds_lon(ij,4), bounds_lat(ij,4))
156            CALL xyz2lonlat(xyz_v(ij+z_down,:), bounds_lon(ij,5), bounds_lat(ij,5))
157            CALL xyz2lonlat(xyz_v(ij+z_rdown,:), bounds_lon(ij,6), bounds_lat(ij,6))
158         ENDDO
159       ENDDO           
160    ENDDO
161         
162!$OMP PARALLEL
163    CALL initialize_physics_omp
164!$OMP END PARALLEL           
165
166    CALL xios_set_context   
167
168
169     
170
171  END SUBROUTINE initialize_physics
172
173
174  SUBROUTINE initialize_physics_omp
175  USE distrib_icosa_lmdz_mod, ONLY : init_distrib_icosa_lmdz, transfer_icosa_to_lmdz
176! from dynamico
177  USE domain_mod
178  USE dimensions
179  USE mpi_mod
180  USE mpipara
181  USE disvert_mod
182  USE xios_mod
183  USE time_mod , ONLY: init_time_icosa=> init_time, dt, itaumax, itau_physics
184  USE omp_para
185
186! from LMDZ
187  USE mod_grid_phy_lmdz, ONLY : unstructured
188  USE mod_phys_lmdz_omp_data, ONLY: klon_omp
189  USE time_phylmdz_mod, ONLY: init_time_lmdz => init_time
190  USE transfert_mod
191  USE physics_distribution_mod, ONLY : init_physics_distribution
192  USE dimphy, ONLY: init_dimphy
193  USE geometry_mod, ONLY : init_geometry
194  USE vertical_layers_mod, ONLY : init_vertical_layers
195  USE planete_mod, ONLY: ini_planete_mod
196!  USE infotrac_phy, ONLY : init_infotrac_phy
197  USE inifis_mod, ONLY : inifis
198!  USE phyaqua_mod, ONLY : iniaqua
199   
200 
201  IMPLICIT NONE
202
203
204
205  INTEGER  :: ind,i,j,ij,pos
206  REAL(rstd),POINTER :: bounds_lon(:,:)
207  REAL(rstd),POINTER :: bounds_lat(:,:)
208 
209  REAL(rstd),ALLOCATABLE :: latfi(:)
210  REAL(rstd),ALLOCATABLE :: lonfi(:)
211  REAL(rstd),ALLOCATABLE :: airefi(:)
212  REAL(rstd),ALLOCATABLE :: bounds_latfi(:,:)
213  REAL(rstd),ALLOCATABLE :: bounds_lonfi(:,:)
214  REAL(rstd),ALLOCATABLE :: ind_cell_glo(:)
215
216  REAL(rstd) :: pseudoalt(llm)
217  REAL(rstd) :: aps(llm)
218  REAL(rstd) :: bps(llm)
219  real(rstd) :: scaleheight
220
221  INTEGER :: run_length 
222  INTEGER :: annee_ref 
223  INTEGER :: day_ref   
224  INTEGER :: day_ini   
225  REAL    :: start_time
226  REAL    :: physics_timestep   
227
228
229!  INTEGER                       :: nqo, nbtr
230  CHARACTER(len=4)              :: type_trac
231!  CHARACTER(len=20),ALLOCATABLE :: tname(:)    ! tracer short name for restart and diagnostics
232!  CHARACTER(len=23),ALLOCATABLE :: ttext(:)     ! tracer long name for diagnostics
233!  INTEGER,ALLOCATABLE           :: niadv(:)    ! equivalent dyn / physique
234!  INTEGER,ALLOCATABLE           :: conv_flg(:) ! conv_flg(it)=0 : convection desactivated for tracer number it
235!  INTEGER,ALLOCATABLE           :: pbl_flg(:)  ! pbl_flg(it)=0  : boundary layer diffusion desactivaded for tracer number it
236!  CHARACTER(len=8),ALLOCATABLE  :: solsym(:)  ! tracer name from inca
237  TYPE(t_field),POINTER,SAVE    :: f_ind_cell_glo(:)
238 
239  INTEGER :: iflag_phys   
240  INTEGER :: nq
241
242    CALL init_distrib_icosa_lmdz
243   
244    ALLOCATE(latfi(klon_omp))
245    ALLOCATE(lonfi(klon_omp))
246    ALLOCATE(airefi(klon_omp))
247    ALLOCATE(bounds_latfi(klon_omp,6))
248    ALLOCATE(bounds_lonfi(klon_omp,6))
249    ALLOCATE(ind_cell_glo(klon_omp))
250
251    CALL transfer_icosa_to_lmdz(geom%lat_i,latfi)
252    CALL transfer_icosa_to_lmdz(geom%lon_i,lonfi)
253    CALL transfer_icosa_to_lmdz(f_bounds_lat,bounds_latfi)
254    CALL transfer_icosa_to_lmdz(f_bounds_lon,bounds_lonfi)
255    CALL transfer_icosa_to_lmdz(geom%Ai,airefi)
256
257    CALL allocate_field(f_ind_cell_glo,field_t,type_real)
258   
259    DO ind=1,ndomain
260      IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master ) CYCLE
261      CALL swap_dimensions(ind)
262      CALL swap_geometry(ind)
263      DO j=jj_begin,jj_end
264        DO i=ii_begin,ii_end
265          ij=(j-1)*iim+i
266          f_ind_cell_glo(ind)%rval2d(ij)=domain(ind)%assign_cell_glo(i,j)
267        ENDDO
268      ENDDO
269    ENDDO
270
271!$OMP BARRIER
272     
273    CALL transfer_icosa_to_lmdz(f_ind_cell_glo,ind_cell_glo)
274    CALL deallocate_field(f_ind_cell_glo)
275     
276             
277    ! Initialize dimphy module
278    CALL init_dimphy(klon_omp,llm)
279
280    CALL init_geometry(klon_omp, lonfi, latfi, bounds_lonfi, bounds_latfi, airefi, INT(ind_cell_glo))
281
282    scaleheight=scale_height/1000. ! Atmospheric scale height (km)
283    aps(1:llm)=0.5*(ap(1:llm)+ap(2:llm+1))
284    bps(1:llm)=0.5*(bp(1:llm)+bp(2:llm+1))
285    pseudoalt(:)=-scaleheight*log(presnivs(:)/preff)
286    CALL init_vertical_layers(llm,preff,scaleheight,ap,bp,aps,bps,presnivs,pseudoalt)
287
288    ! Initialize planet_mod (quite redundant wrt vertical_levels...)
289    CALL ini_planete_mod(llm,preff,ap,bp)
290
291    ! Initialize tracer names, numbers, etc. for physics
292!$OMP MASTER
293    !Config  Key  = type_trac
294    !Config  Desc = Choix de couplage avec model de chimie INCA ou REPROBUS
295    !Config  Def  = lmdz
296    !Config  Help =
297    !Config         'lmdz' = pas de couplage, pur LMDZ
298    !Config         'inca' = model de chime INCA
299    !Config         'repr' = model de chime REPROBUS
300     type_trac = 'lmdz'
301     CALL getin('type_trac',type_trac)
302
303! init model for standard lmdz case
304!    nqo=2
305!    nbtr=2
306    ALLOCATE(tname(nqtot))
307!    ALLOCATE(ttext(nqtot))
308!    ALLOCATE(niadv(nqtot))
309!    ALLOCATE(conv_flg(nbtr))
310!    ALLOCATE(pbl_flg(nbtr))
311!    ALLOCATE(solsym(nbtr))
312   
313
314! read tname() from traceur.def file
315    OPEN(unit=42,file="traceur.def",form="formatted",status="old",iostat=ierr)
316    IF (ierr==0) THEN
317      READ(42,*) nq ! should be the same as nqtot
318      IF (nq /= nqtot) THEN
319        WRITE(*,*) "Error: number of tracers in tracer.def should match nqtot!"
320        WRITE(*,*) "       will just use nqtot=",nqtot," tracers"
321      ENDIF
322      DO i=1,nqtot
323        READ(42,*) tname(i)
324      ENDDO
325      CLOSE(42)
326    ENDIF
327!$OMP END MASTER
328!$OMP BARRIER
329
330!    CALL init_infotrac_phy(nqtot,nqo,nbtr,tname,ttext,type_trac,&
331!                           niadv,conv_flg,pbl_flg,solsym)
332
333   ! Initialize physical constant
334!    day_length=86400
335!    CALL getin('day_length',day_length)
336!    CALL inifis(day_length,radius,g,kappa*cpp,cpp)
337
338   ! Initialize physical constant
339!    CALL inifis(klon_omp,nlayer,nqtot,pdayref,punjours,nday,ptimestep, &
340!            latitude,longitude,cell_area,prad,pg,pr,pcpp)
341    start_day=0
342    CALL getin('start_day',start_day)
343    day_length=86400
344    CALL getin('day_length',day_length)
345    ndays=nint(itaumax*(dt/day_length))! number of days to run
346    physics_timestep=dt*itau_physics
347    CALL inifis(klon_omp,llm,nqtot,start_day,day_length,ndays,physics_timestep, &
348            latfi,lonfi,airefi,radius,g,kappa*cpp,cpp)
349
350   
351  ! init time
352!    annee_ref=2015
353!    CALL getin("anneeref",annee_ref)
354   
355!    day_ref=1
356!    CALL getin("dayref",day_ref)
357   
358!    physics_timestep=dt*itau_physics
359!    run_length=itaumax*dt
360!    ndays=NINT(run_length/day_length)
361   
362!    day_ini=INT(itau0*dt/day_length)+day_ref
363!    start_time= itau0*dt/day_length-INT(itau0*dt/day_length)
364
365!    CALL init_time_lmdz(annee_ref, day_ref, day_ini, start_time, ndays, physics_timestep)
366
367!  Additional initializations for aquaplanets
368!    CALL getin("iflag_phys",iflag_phys)
369!    IF (iflag_phys>=100) THEN
370!      CALL iniaqua(klon_omp, iflag_phys)
371!    END IF
372
373    ! Initialize "calendar"
374!$OMP MASTER
375    day_length=86400
376    CALL getin('day_length',day_length)
377    year_length = 31557600 ! 365.25 * 86400
378    CALL getin('year_length',year_length)
379    start_day = 0
380    CALL getin('start_day',start_day)
381    ! initialize pday and ptime
382    pday = start_day
383    ptime = 0 ! assume we start at midnight ...
384!$OMP END MASTER
385!$OMP BARRIER
386 
387  END SUBROUTINE  initialize_physics_omp
388 
389 
390
391
392  SUBROUTINE physics
393  USE icosa
394  USE time_mod
395  USE disvert_mod
396  USE transfert_mod
397  USE mpipara
398  USE xios_mod
399  USE wxios
400  USE trace
401  USE distrib_icosa_lmdz_mod, ONLY : transfer_icosa_to_lmdz, transfer_lmdz_to_icosa
402  USE physics_external_mod, ONLY : it, f_phis, f_ps, f_theta_rhodz, f_u, f_wflux, f_q
403  USE write_field_mod
404  USE checksum_mod
405! from LMDZ
406  USE mod_phys_lmdz_omp_data, ONLY: klon_omp
407  USE geometry_mod, ONLY : cell_area
408  USE physiq_mod, ONLY: physiq
409  IMPLICIT NONE
410 
411    REAL(rstd),POINTER :: phis(:)
412    REAL(rstd),POINTER :: ps(:)
413    REAL(rstd),POINTER :: theta_rhodz(:,:,:)
414    REAL(rstd),POINTER :: u(:,:)
415    REAL(rstd),POINTER :: wflux(:,:)
416    REAL(rstd),POINTER :: q(:,:,:)
417    REAL(rstd),POINTER :: p(:,:)
418    REAL(rstd),POINTER :: pks(:)
419    REAL(rstd),POINTER :: pk(:,:)
420    REAL(rstd),POINTER :: p_layer(:,:)
421    REAL(rstd),POINTER :: theta(:,:)
422    REAL(rstd),POINTER :: phi(:,:)
423    REAL(rstd),POINTER :: Temp(:,:)
424    REAL(rstd),POINTER :: ulon(:,:)
425    REAL(rstd),POINTER :: ulat(:,:)
426    REAL(rstd),POINTER :: dulon(:,:)
427    REAL(rstd),POINTER :: dulat(:,:)
428    REAL(rstd),POINTER :: dTemp(:,:)
429    REAL(rstd),POINTER :: dq(:,:,:)
430    REAL(rstd),POINTER :: dps(:)
431    REAL(rstd),POINTER :: duc(:,:,:)
432
433
434    INTEGER :: ind,l
435   
436    REAL(rstd),ALLOCATABLE,SAVE :: ps_phy(:)
437!$OMP THREADPRIVATE(ps_phy)
438    REAL(rstd),ALLOCATABLE,SAVE :: p_phy(:,:)
439!$OMP THREADPRIVATE(p_phy)
440    REAL(rstd),ALLOCATABLE,SAVE :: p_layer_phy(:,:)
441!$OMP THREADPRIVATE(p_layer_phy)
442    REAL(rstd),ALLOCATABLE,SAVE :: Temp_phy(:,:)
443!$OMP THREADPRIVATE(Temp_phy)
444    REAL(rstd),ALLOCATABLE,SAVE :: phis_phy(:)
445!$OMP THREADPRIVATE(phis_phy)
446    REAL(rstd),ALLOCATABLE,SAVE :: phi_phy(:,:)
447!$OMP THREADPRIVATE(phi_phy)
448    REAL(rstd),ALLOCATABLE,SAVE :: ulon_phy(:,:)
449!$OMP THREADPRIVATE(ulon_phy)
450    REAL(rstd),ALLOCATABLE,SAVE :: ulat_phy(:,:)
451!$OMP THREADPRIVATE(ulat_phy)
452    REAL(rstd),ALLOCATABLE,SAVE :: q_phy(:,:,:)
453!$OMP THREADPRIVATE(q_phy)
454    REAL(rstd),ALLOCATABLE,SAVE :: wflux_phy(:,:)
455!$OMP THREADPRIVATE(wflux_phy)
456    REAL(rstd),ALLOCATABLE,SAVE :: dulon_phy(:,:)
457!$OMP THREADPRIVATE(dulon_phy)
458    REAL(rstd),ALLOCATABLE,SAVE :: dulat_phy(:,:)
459!$OMP THREADPRIVATE(dulat_phy)
460    REAL(rstd),ALLOCATABLE,SAVE :: dTemp_phy(:,:)
461!$OMP THREADPRIVATE(dTemp_phy)
462    REAL(rstd),ALLOCATABLE,SAVE :: dq_phy(:,:,:)
463!$OMP THREADPRIVATE(dq_phy)
464    REAL(rstd),ALLOCATABLE,SAVE :: dps_phy(:)
465!$OMP THREADPRIVATE(dps_phy)
466    REAL(rstd)   :: dtphy
467    LOGICAL      :: debut
468    LOGICAL      :: lafin
469    LOGICAL,SAVE :: first=.TRUE.
470!$OMP THREADPRIVATE(first)
471
472   
473    IF(first) THEN
474      debut=.TRUE.
475    ELSE
476      debut=.FALSE.
477    ENDIF
478
479
480    IF(it-itau0>=itaumax) THEN
481      lafin=.TRUE.
482    ELSE
483      lafin=.FALSE.
484    ENDIF
485
486    IF (first) THEN
487      first=.FALSE.
488      CALL init_message(f_u,req_e1_vect,req_u)
489      ALLOCATE(ps_phy(klon_omp))
490      ALLOCATE(p_phy(klon_omp,llm+1))
491      ALLOCATE(p_layer_phy(klon_omp,llm))
492      ALLOCATE(Temp_phy(klon_omp,llm))
493      ALLOCATE(phis_phy(klon_omp))
494      ALLOCATE(phi_phy(klon_omp,llm))
495      ALLOCATE(ulon_phy(klon_omp,llm))
496      ALLOCATE(ulat_phy(klon_omp,llm))
497      ALLOCATE(q_phy(klon_omp,llm,nqtot))
498      ALLOCATE(wflux_phy(klon_omp,llm))
499      ALLOCATE(dulon_phy(klon_omp,llm))
500      ALLOCATE(dulat_phy(klon_omp,llm))
501      ALLOCATE(dTemp_phy(klon_omp,llm))
502      ALLOCATE(dq_phy(klon_omp,llm,nqtot))
503      ALLOCATE(dps_phy(klon_omp))
504!$OMP BARRIER
505    ENDIF
506
507
508!$OMP MASTER       
509!    CALL update_calendar(it)
510!$OMP END MASTER
511!$OMP BARRIER
512    dtphy=itau_physics*dt
513   
514   
515   
516    CALL transfert_message(f_u,req_u)
517   
518    DO ind=1,ndomain
519      CALL swap_dimensions(ind)
520      IF (assigned_domain(ind)) THEN
521        CALL swap_geometry(ind)
522     
523        phis=f_phis(ind)
524        ps=f_ps(ind)
525        theta_rhodz=f_theta_rhodz(ind)
526        u=f_u(ind)
527        q=f_q(ind)
528        wflux=f_wflux(ind)
529        p=f_p(ind)
530        pks=f_pks(ind)
531        pk=f_pk(ind)
532        p_layer=f_p_layer(ind)
533        theta=f_theta(ind)
534        phi=f_phi(ind)
535        Temp=f_Temp(ind)
536        ulon=f_ulon(ind)
537        ulat=f_ulat(ind)
538           
539        CALL grid_icosa_to_physics
540
541      ENDIF
542    ENDDO
543   
544!$OMP BARRIER
545!$OMP MASTER
546    CALL SYSTEM_CLOCK(start_clock)
547!$OMP END MASTER
548    CALL trace_start("physic")
549!    CALL trace_off()
550
551
552!    CALL writeField("p_in",f_p)
553!    CALL writeField("p_layer_in",f_p_layer)
554!    CALL writeField("phi_in",f_phi)
555!    CALL writeField("phis_in",f_phis)
556!    CALL writeField("ulon_in",f_ulon)
557!    CALL writeField("ulat_in",f_ulat)
558!    CALL writeField("Temp_in",f_Temp)
559!    CALL writeField("q_in",f_q)
560!    CALL writeField("wflux_in",f_wflux)
561
562!    CALL checksum(f_p)
563!    CALL checksum(f_p_layer)
564!    CALL checksum(f_phi)
565!    CALL checksum(f_phis)
566!    CALL checksum(f_ulon)
567!    CALL checksum(f_ulat)
568!    CALL checksum(f_Temp)
569!    CALL checksum(f_q)
570!    CALL checksum(f_wflux)
571
572    CALL transfer_icosa_to_lmdz(f_p      , p_phy)
573    CALL transfer_icosa_to_lmdz(f_p_layer, p_layer_phy)
574    CALL transfer_icosa_to_lmdz(f_phi    , phi_phy)
575    CALL transfer_icosa_to_lmdz(f_phis   , phis_phy )
576    CALL transfer_icosa_to_lmdz(f_ulon   , ulon_phy )
577    CALL transfer_icosa_to_lmdz(f_ulat   , ulat_phy)
578    CALL transfer_icosa_to_lmdz(f_Temp   , Temp_phy)
579    CALL transfer_icosa_to_lmdz(f_q      , q_phy)
580    CALL transfer_icosa_to_lmdz(f_wflux  , wflux_phy)
581
582    DO l=1,llm
583      wflux_phy(:,l)=wflux_phy(:,l)*cell_area(:)
584      phi_phy(:,l)=phi_phy(:,l)-phis_phy(:)
585    ENDDO
586   
587    CALL wxios_set_context()
588
589    ! Update pday and ptime to send to physics
590!$OMP MASTER
591    ptime=ptime+dtphy/day_length
592    IF (ptime >= 1.) THEN
593      ptime=ptime-1
594      pday=pday+1
595    ENDIF
596!$OMP END MASTER
597!$OMP BARRIER   
598    CALL physiq(klon_omp, llm, nqtot, &
599                tname, &
600                debut, lafin, &
601                pday, ptime, dtphy, &
602                p_phy, p_layer_phy, phi_phy, &
603                ulon_phy, ulat_phy, Temp_phy, q_phy, &
604                wflux_phy, &
605                dulon_phy, dulat_phy, dTemp_phy, dq_phy, dps_phy)
606   
607    CALL transfer_lmdz_to_icosa(dulon_phy, f_dulon )
608    CALL transfer_lmdz_to_icosa(dulat_phy, f_dulat )
609    CALL transfer_lmdz_to_icosa(dTemp_phy, f_dTemp )
610    CALL transfer_lmdz_to_icosa(dq_phy   , f_dq )
611    CALL transfer_lmdz_to_icosa(dps_phy  , f_dps )
612 
613!    CALL writeField("dulon_out",f_dulon)
614!    CALL writeField("dulat_out",f_dulat)
615!    CALL writeField("dTemp_out",f_dTemp)
616!    CALL writeField("dq_out",f_dq)
617!    CALL writeField("dps_out",f_dps)
618
619!    CALL checksum(f_dulon)
620!    CALL checksum(f_dulat)
621!    CALL checksum(f_dTemp)
622!    CALL checksum(f_dq)
623!    CALL checksum(f_dps)
624   
625    CALL send_message(f_dps,req_dps0)
626    CALL send_message(f_dulon,req_dulon0)
627    CALL send_message(f_dulat,req_dulat0)
628    CALL send_message(f_dTemp,req_dTemp0)
629    CALL send_message(f_dq,req_dq0)
630
631    CALL wait_message(req_dps0)
632    CALL wait_message(req_dulon0)
633    CALL wait_message(req_dulat0)
634    CALL wait_message(req_dTemp0)
635    CALL wait_message(req_dq0)
636
637
638!    CALL trace_on()
639    CALL trace_end("physic")
640!$OMP MASTER
641    CALL SYSTEM_CLOCK(stop_clock)
642    count_clock=count_clock+stop_clock-start_clock
643!$OMP END MASTER
644
645!$OMP BARRIER                       
646
647    DO ind=1,ndomain
648      CALL swap_dimensions(ind)
649      IF (assigned_domain(ind)) THEN
650        CALL swap_geometry(ind)
651
652        theta_rhodz=f_theta_rhodz(ind)
653        u=f_u(ind)
654        q=f_q(ind)
655        ps=f_ps(ind)
656        dulon=f_dulon(ind)
657        dulat=f_dulat(ind)
658        Temp=f_temp(ind)
659        dTemp=f_dTemp(ind)
660        dq=f_dq(ind)
661        dps=f_dps(ind)
662        duc=f_duc(ind)
663        p=f_p(ind)
664        pks=f_pks(ind)
665        pk=f_pk(ind)
666     
667        CALL grid_physics_to_icosa
668      ENDIF
669    ENDDO
670
671!$OMP BARRIER
672    CALL xios_set_context   
673   
674 
675  CONTAINS
676
677    SUBROUTINE grid_icosa_to_physics
678    USE pression_mod
679    USE exner_mod
680    USE theta2theta_rhodz_mod
681    USE geopotential_mod
682    USE wind_mod
683    USE omp_para
684    IMPLICIT NONE
685   
686    REAL(rstd) :: uc(3)
687    INTEGER :: i,j,ij,l
688   
689
690! compute pression
691
692      DO    l    = ll_begin,ll_endp1
693        DO j=jj_begin,jj_end
694          DO i=ii_begin,ii_end
695            ij=(j-1)*iim+i
696            p(ij,l) = ap(l) + bp(l) * ps(ij)
697          ENDDO
698        ENDDO
699      ENDDO
700
701!$OMP BARRIER
702
703! compute exner
704       
705       IF (is_omp_first_level) THEN
706         DO j=jj_begin,jj_end
707            DO i=ii_begin,ii_end
708               ij=(j-1)*iim+i
709               pks(ij) = cpp * ( ps(ij)/preff ) ** kappa
710            ENDDO
711         ENDDO
712       ENDIF
713
714       ! 3D : pk
715       DO l = ll_begin,ll_end
716          DO j=jj_begin,jj_end
717             DO i=ii_begin,ii_end
718                ij=(j-1)*iim+i
719                pk(ij,l) = cpp * ((.5/preff)*(p(ij,l)+p(ij,l+1))) ** kappa
720             ENDDO
721          ENDDO
722       ENDDO
723
724!$OMP BARRIER
725
726!   compute theta, temperature and pression at layer
727    DO    l    = ll_begin, ll_end
728      DO j=jj_begin,jj_end
729        DO i=ii_begin,ii_end
730          ij=(j-1)*iim+i
731          theta(ij,l) = theta_rhodz(ij,l,1) / ((p(ij,l)-p(ij,l+1))/g)
732          Temp(ij,l) = theta(ij,l) * pk(ij,l) / cpp
733          p_layer(ij,l)=preff*(pk(ij,l)/cpp)**(1./kappa)
734        ENDDO
735      ENDDO
736    ENDDO
737
738
739!!! Compute geopotential
740       
741  ! for first layer
742  IF (is_omp_first_level) THEN
743    DO j=jj_begin,jj_end
744      DO i=ii_begin,ii_end
745        ij=(j-1)*iim+i
746        phi( ij,1 ) = phis( ij ) + theta(ij,1) * ( pks(ij) - pk(ij,1) )
747      ENDDO
748    ENDDO
749  ENDIF
750!!-> implicit flush on phi(:,1)
751         
752  ! for other layers
753  DO l = ll_beginp1, ll_end
754    DO j=jj_begin,jj_end
755      DO i=ii_begin,ii_end
756        ij=(j-1)*iim+i
757        phi(ij,l) =  0.5 * ( theta(ij,l)  + theta(ij,l-1) )  &
758                         * (  pk(ij,l-1) -  pk(ij,l)    )
759      ENDDO
760    ENDDO
761  ENDDO       
762
763!$OMP BARRIER
764
765
766  IF (is_omp_first_level) THEN
767    DO l = 2, llm
768      DO j=jj_begin,jj_end
769! ---> Bug compilo intel ici en openmp
770! ---> Couper la boucle
771       IF (j==jj_end+1) PRINT*,"this message must not be printed"
772        DO i=ii_begin,ii_end
773          ij=(j-1)*iim+i
774          phi(ij,l) = phi(ij,l)+ phi(ij,l-1)
775        ENDDO
776      ENDDO
777    ENDDO
778! --> IMPLICIT FLUSH on phi --> non
779  ENDIF
780
781! compute wind centered lon lat compound
782    DO l=ll_begin,ll_end
783      DO j=jj_begin,jj_end
784        DO i=ii_begin,ii_end
785          ij=(j-1)*iim+i
786          uc(:)=1/Ai(ij)*                                                                                                &
787                        ( ne(ij,right)*u(ij+u_right,l)*le(ij+u_right)*((xyz_v(ij+z_rdown,:)+xyz_v(ij+z_rup,:))/2-centroid(ij,:))  &
788                         + ne(ij,rup)*u(ij+u_rup,l)*le(ij+u_rup)*((xyz_v(ij+z_rup,:)+xyz_v(ij+z_up,:))/2-centroid(ij,:))          &
789                         + ne(ij,lup)*u(ij+u_lup,l)*le(ij+u_lup)*((xyz_v(ij+z_up,:)+xyz_v(ij+z_lup,:))/2-centroid(ij,:))          &
790                         + ne(ij,left)*u(ij+u_left,l)*le(ij+u_left)*((xyz_v(ij+z_lup,:)+xyz_v(ij+z_ldown,:))/2-centroid(ij,:))    &
791                         + ne(ij,ldown)*u(ij+u_ldown,l)*le(ij+u_ldown)*((xyz_v(ij+z_ldown,:)+xyz_v(ij+z_down,:))/2-centroid(ij,:))&
792                         + ne(ij,rdown)*u(ij+u_rdown,l)*le(ij+u_rdown)*((xyz_v(ij+z_down,:)+xyz_v(ij+z_rdown,:))/2-centroid(ij,:)))
793          ulon(ij,l)=sum(uc(:)*elon_i(ij,:))
794          ulat(ij,l)=sum(uc(:)*elat_i(ij,:))
795        ENDDO
796      ENDDO
797    ENDDO
798
799!$OMP BARRIER
800    END SUBROUTINE grid_icosa_to_physics
801
802
803    SUBROUTINE grid_physics_to_icosa
804    USE theta2theta_rhodz_mod
805    USE omp_para
806    IMPLICIT NONE
807      INTEGER :: i,j,ij,l,iq
808         
809      DO l=ll_begin,ll_end
810        DO j=jj_begin,jj_end
811          DO i=ii_begin,ii_end
812            ij=(j-1)*iim+i
813            duc(ij,:,l)=dulon(ij,l)*elon_i(ij,:)+dulat(ij,l)*elat_i(ij,:)
814          ENDDO
815        ENDDO
816      ENDDO
817
818      DO l=ll_begin,ll_end
819        DO j=jj_begin,jj_end
820          DO i=ii_begin,ii_end
821            ij=(j-1)*iim+i
822            u(ij+u_right,l) = u(ij+u_right,l) + dtphy * sum( 0.5*(duc(ij,:,l) + duc(ij+t_right,:,l))*ep_e(ij+u_right,:) )
823            u(ij+u_lup,l) = u(ij+u_lup,l) + dtphy * sum( 0.5*(duc(ij,:,l) + duc(ij+t_lup,:,l))*ep_e(ij+u_lup,:) )
824            u(ij+u_ldown,l) = u(ij+u_ldown,l) + dtphy*sum( 0.5*(duc(ij,:,l) + duc(ij+t_ldown,:,l))*ep_e(ij+u_ldown,:) )
825          ENDDO
826        ENDDO
827      ENDDO         
828
829      DO l=ll_begin,ll_end
830        DO j=jj_begin,jj_end
831          DO i=ii_begin,ii_end
832            ij=(j-1)*iim+i
833            Temp(ij,l)=Temp(ij,l)+ dtphy * dTemp(ij,l)
834          ENDDO
835        ENDDO
836      ENDDO         
837     
838      DO iq=1,nqtot
839        DO l=ll_begin,ll_end
840          DO j=jj_begin,jj_end
841            DO i=ii_begin,ii_end
842              ij=(j-1)*iim+i
843              q(ij,l,iq)=q(ij,l,iq)+ dtphy * dq(ij,l,iq)
844            ENDDO
845          ENDDO
846        ENDDO
847      ENDDO
848
849!$OMP BARRIER
850     
851     IF (is_omp_first_level) THEN
852       DO j=jj_begin,jj_end
853         DO i=ii_begin,ii_end
854           ij=(j-1)*iim+i
855           ps(ij)=ps(ij)+ dtphy * dps(ij)
856          ENDDO
857       ENDDO
858     ENDIF
859
860!     CALL compute_temperature2theta_rhodz(ps,Temp,theta_rhodz,0)
861
862! compute pression
863!$OMP BARRIER
864      DO    l    = ll_begin,ll_endp1
865        DO j=jj_begin,jj_end
866          DO i=ii_begin,ii_end
867            ij=(j-1)*iim+i
868            p(ij,l) = ap(l) + bp(l) * ps(ij)
869          ENDDO
870        ENDDO
871      ENDDO
872
873!$OMP BARRIER
874
875! compute exner
876       
877       IF (is_omp_first_level) THEN
878         DO j=jj_begin,jj_end
879            DO i=ii_begin,ii_end
880               ij=(j-1)*iim+i
881               pks(ij) = cpp * ( ps(ij)/preff ) ** kappa
882            ENDDO
883         ENDDO
884       ENDIF
885
886       ! 3D : pk
887       DO l = ll_begin,ll_end
888          DO j=jj_begin,jj_end
889             DO i=ii_begin,ii_end
890                ij=(j-1)*iim+i
891                pk(ij,l) = cpp * ((.5/preff)*(p(ij,l)+p(ij,l+1))) ** kappa
892             ENDDO
893          ENDDO
894       ENDDO
895
896!$OMP BARRIER
897
898!   compute theta, temperature and pression at layer
899    DO    l    = ll_begin, ll_end
900      DO j=jj_begin,jj_end
901        DO i=ii_begin,ii_end
902          ij=(j-1)*iim+i
903          theta_rhodz(ij,l,1) = temp(ij,l) * ((p(ij,l)-p(ij,l+1))/g) / (pk(ij,l) / cpp )
904        ENDDO
905      ENDDO
906    ENDDO
907   
908    END SUBROUTINE grid_physics_to_icosa
909
910
911
912  END SUBROUTINE physics
913
914
915
916
917
918END MODULE interface_icosa_lmdz_mod
Note: See TracBrowser for help on using the repository browser.