source: ICOSA_LMDZ/src/phylmd/interface_icosa_lmdz.F90 @ 4857

Last change on this file since 4857 was 4857, checked in by yann meurdesoif, 2 months ago

Fix incomptability with recent version of dynamico due to convergence hexa<->unstruct

YM

File size: 31.4 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, req_z
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_vort(:)   
20  TYPE(t_field),POINTER,SAVE :: f_vortc(:)   
21  TYPE(t_field),POINTER,SAVE :: f_dulon(:)
22  TYPE(t_field),POINTER,SAVE :: f_dulat(:)
23  TYPE(t_field),POINTER,SAVE :: f_dTemp(:)
24  TYPE(t_field),POINTER,SAVE :: f_dq(:)
25  TYPE(t_field),POINTER,SAVE :: f_dps(:)
26  TYPE(t_field),POINTER,SAVE :: f_duc(:)
27  TYPE(t_field),POINTER,SAVE :: f_bounds_lon(:)
28  TYPE(t_field),POINTER,SAVE :: f_bounds_lat(:)
29
30  INTEGER :: start_clock
31  INTEGER :: stop_clock
32  INTEGER :: count_clock=0
33 
34  INTEGER,SAVE :: nbp_phys
35  INTEGER,SAVE :: nbp_phys_glo
36
37
38CONTAINS
39
40  SUBROUTINE pre_initialize_physics
41  USE etat0_plugin_mod
42  USE icolmdz_etat0
43  USE icolmdz_param_gravity_wave, ONLY: init_param_gravity_wave
44  IMPLICIT NONE
45  !$OMP PARALLEL
46    init_etat0_plugin => init_etat0
47    etat0_plugin => etat0
48    CALL init_param_gravity_wave
49  !$OMP END PARALLEL
50  END SUBROUTINE pre_initialize_physics
51
52  SUBROUTINE initialize_physics
53  USE distrib_icosa_lmdz_mod, ONLY : init_distrib_icosa_lmdz, transfer_icosa_to_lmdz
54! from dynamico
55  USE domain_mod
56  USE dimensions
57  USE mpi_mod
58  USE mpipara
59  USE disvert_mod
60  USE xios_mod
61  USE time_mod , init_time_icosa=> init_time
62  USE transfert_mod
63  USE nudging_mod, ONLY : lam_halo_scheme
64 
65! from LMDZ
66  USE mod_grid_phy_lmdz, ONLY : unstructured
67  USE mod_phys_lmdz_omp_data, ONLY: klon_omp
68  USE transfert_mod
69  USE physics_distribution_mod, ONLY : init_physics_distribution
70  USE infotrac_phy, ONLY: init_infotrac_phy
71  USE icolmdz_param_gravity_wave, ONLY: param_gravity_wave
72
73   
74 
75  IMPLICIT NONE
76  INTEGER  :: ind,i,j,ij,pos,h
77  REAL(rstd),POINTER :: bounds_lon(:,:)
78  REAL(rstd),POINTER :: bounds_lat(:,:)
79 
80  REAL(rstd),ALLOCATABLE :: latfi(:)
81  REAL(rstd),ALLOCATABLE :: lonfi(:)
82  REAL(rstd),ALLOCATABLE :: airefi(:)
83  REAL(rstd),ALLOCATABLE :: bounds_latfi(:,:)
84  REAL(rstd),ALLOCATABLE :: bounds_lonfi(:,:)
85  LOGICAL   ,ALLOCATABLE :: outside(:,:)
86  LOGICAL   ,ALLOCATABLE :: outside_tmp(:,:)
87  LOGICAL   ,POINTER     :: out(:,:)
88!  REAL(rstd) :: pseudoalt(llm)
89
90  INTEGER :: nbp_phys, nbp_phys_glo
91 
92!$OMP PARALLEL
93    CALL allocate_field(f_bounds_lon,field_t,type_real,6)
94    CALL allocate_field(f_bounds_lat,field_t,type_real,6)
95    CALL allocate_field(f_p,field_t,type_real,llm+1,name="p_in")
96    CALL allocate_field(f_pks,field_t,type_real)
97    CALL allocate_field(f_pk,field_t,type_real,llm)
98    CALL allocate_field(f_p_layer,field_t,type_real,llm,name="p_layer_in")
99    CALL allocate_field(f_theta,field_t,type_real,llm)
100    CALL allocate_field(f_phi,field_t,type_real,llm,name="phi_in")
101    CALL allocate_field(f_Temp,field_t,type_real,llm,name="Temp_in")
102    CALL allocate_field(f_ulon,field_t,type_real,llm,name="ulon_in")
103    CALL allocate_field(f_ulat,field_t,type_real,llm,name="ulat_in")
104    CALL allocate_field(f_vort,field_z,type_real,llm,name="vort_in")
105    CALL allocate_field(f_vortc,field_t,type_real,llm,name="vortc_in")
106    CALL allocate_field(f_dulon,field_t,type_real,llm,name="dulon_out")
107    CALL allocate_field(f_dulat,field_t,type_real,llm,name="dulat_out")
108    CALL allocate_field(f_dTemp,field_t,type_real,llm,name="dTemp_out")
109    CALL allocate_field(f_dq,field_t,type_real,llm,nqtot,name="dq_out")
110    CALL allocate_field(f_dps,field_t,type_real,name="dps_out")
111    CALL allocate_field(f_duc,field_t,type_real,3,llm)   
112
113    CALL init_message(f_dps,req_i0,req_dps0)
114    CALL init_message(f_dulon,req_i0,req_dulon0)
115    CALL init_message(f_dulat,req_i0,req_dulat0)
116    CALL init_message(f_dTemp,req_i0,req_dTemp0)
117    CALL init_message(f_dq,req_i0,req_dq0)
118!$OMP END PARALLEL   
119
120    nbp_phys=0
121    DO ind=1,ndomain
122      CALL swap_dimensions(ind)
123
124      ALLOCATE(outside(ii_begin:ii_end,jj_begin:jj_end)) ! for limited area : don't take cells arround the border
125      ALLOCATE(outside_tmp(ii_begin-1:ii_end+1,jj_begin-1:jj_end+1)) ! for limited area : don't take cells arround the border
126      out=>domain(ind)%outside
127      DO j=jj_begin,jj_end
128        DO i=ii_begin,ii_end
129          outside(i,j)=  out(i+1,j)     .OR. & ! right
130                         out(i,j+1    ) .OR. & ! rup
131                         out(i-1  ,j+1) .OR. & ! lup
132                         out(i-1  ,j)   .OR. & !left
133                         out(i    ,j-1) .OR. & !ldown
134                         out(i+1,j-1)          !rdown   
135        ENDDO
136      ENDDO
137
138      outside_tmp=.FALSE.
139      outside_tmp(ii_begin:ii_end,jj_begin:jj_end)=outside
140     
141      DO h=1,lam_halo_scheme-1 ! do not compute physic on limited area halo
142        DO j=jj_begin,jj_end
143          DO i=ii_begin,ii_end
144              outside(i,j) = outside_tmp(i,j)       .OR. &
145                             outside_tmp(i+1,j)     .OR. & ! right
146                             outside_tmp(i,j+1    ) .OR. & ! rup
147                             outside_tmp(i-1  ,j+1) .OR. & ! lup
148                             outside_tmp(i-1  ,j)   .OR. & !left
149                             outside_tmp(i    ,j-1) .OR. & !ldown
150                             outside_tmp(i+1,j-1)          !rdown
151           ENDDO
152        ENDDO
153        outside_tmp(ii_begin:ii_end,jj_begin:jj_end)=outside
154      ENDDO
155     
156      DO j=jj_begin,jj_end
157        DO i=ii_begin,ii_end
158          IF (domain(ind)%own(i,j) .AND. .NOT.outside(i,j)) nbp_phys=nbp_phys+1
159        ENDDO
160      ENDDO
161      DEALLOCATE(outside)
162      DEALLOCATE(outside_tmp)
163    ENDDO
164   
165
166!initialize LMDZ5 physic mpi decomposition
167    CALL MPI_ALLREDUCE(nbp_phys,nbp_phys_glo,1,MPI_INTEGER,MPI_SUM,comm_icosa,ierr)
168    CALL init_physics_distribution(unstructured, 6, nbp_phys, 1, nbp_phys_glo, llm, comm_icosa)
169   
170    DO ind=1,ndomain
171        CALL swap_dimensions(ind)
172        CALL swap_geometry(ind)
173        bounds_lon=f_bounds_lon(ind)
174        bounds_lat=f_bounds_lat(ind)
175        DO j=jj_begin,jj_end
176          DO i=ii_begin,ii_end
177            ij=(j-1)*iim+i
178            CALL xyz2lonlat(xyz_v(ij+z_rup,:), bounds_lon(ij,1), bounds_lat(ij,1))
179            CALL xyz2lonlat(xyz_v(ij+z_up,:), bounds_lon(ij,2), bounds_lat(ij,2))
180            CALL xyz2lonlat(xyz_v(ij+z_lup,:), bounds_lon(ij,3), bounds_lat(ij,3))
181            CALL xyz2lonlat(xyz_v(ij+z_ldown,:), bounds_lon(ij,4), bounds_lat(ij,4))
182            CALL xyz2lonlat(xyz_v(ij+z_down,:), bounds_lon(ij,5), bounds_lat(ij,5))
183            CALL xyz2lonlat(xyz_v(ij+z_rdown,:), bounds_lon(ij,6), bounds_lat(ij,6))
184         ENDDO
185       ENDDO           
186    ENDDO
187
188  CALL init_infotrac_phy
189
190         
191!$OMP PARALLEL
192    CALL initialize_physics_omp
193    CALL param_gravity_wave
194!$OMP END PARALLEL           
195
196    CALL xios_set_context   
197
198  END SUBROUTINE initialize_physics
199
200
201  SUBROUTINE initialize_physics_omp
202  USE distrib_icosa_lmdz_mod, ONLY : init_distrib_icosa_lmdz, transfer_icosa_to_lmdz
203! from dynamico
204  USE domain_mod
205  USE dimensions
206  USE mpi_mod
207  USE mpipara
208  USE disvert_mod
209  USE earth_const, ONLY: scale_height
210  USE xios_mod
211  USE time_mod , init_time_icosa=> init_time
212  USE omp_para
213
214! from LMDZ
215  USE mod_grid_phy_lmdz, ONLY : unstructured, klon_glo
216  USE mod_phys_lmdz_para, ONLY: klon_omp, reduce_min_lmdz => reduce_min , gather_lmdz => gather , bcast_lmdz => bcast
217  USE time_phylmdz_mod, ONLY: init_time_lmdz => init_time
218  USE transfert_mod
219!  USE physics_distribution_mod, ONLY : init_physics_distribution
220  USE geometry_mod, ONLY : init_geometry
221  USE vertical_layers_mod, ONLY : init_vertical_layers
222  USE infotrac_phy, ONLY : init_infotrac_phy
223  USE inifis_mod, ONLY : inifis
224  USE readTracFiles_mod, ONLY: trac_type, isot_type
225  USE tracer_icosa_mod, ONLY : tracs
226   USE readTracFiles_mod, ONLY: delPhase
227!  USE phyaqua_mod, ONLY : iniaqua
228   
229 
230  IMPLICIT NONE
231
232
233
234  INTEGER  :: ind,i,j,k,ij,pos
235  REAL(rstd),POINTER :: bounds_lon(:,:)
236  REAL(rstd),POINTER :: bounds_lat(:,:)
237 
238  REAL(rstd),ALLOCATABLE :: latfi(:)
239  REAL(rstd),ALLOCATABLE :: lonfi(:)
240  REAL(rstd),ALLOCATABLE :: airefi(:)
241  REAL(rstd),ALLOCATABLE :: bounds_latfi(:,:)
242  REAL(rstd),ALLOCATABLE :: bounds_lonfi(:,:)
243  REAL(rstd),ALLOCATABLE :: ind_cell_glo_r(:)
244  INTEGER,   ALLOCATABLE :: ind_cell_glo(:)
245  INTEGER,   ALLOCATABLE :: ind_cell_glo_tot(:)
246  INTEGER,   ALLOCATABLE :: cell_glo_tot(:)
247  INTEGER :: ncell_glo_tot
248
249  REAL(rstd) :: pseudoalt(llm)
250  REAL(rstd) :: aps(llm)
251  REAL(rstd) :: bps(llm)
252  REAL(rstd) :: scaleheight
253
254  INTEGER :: run_length 
255  REAL :: day_length ! length of a day (s) ! SAVEd to be OpenMP shared <--- NO!!!!
256  INTEGER :: annee_ref 
257  INTEGER :: day_ref   
258  INTEGER :: day_ini   
259  REAL    :: start_time
260  REAL    :: physics_timestep   
261
262  ! Tracer stuff (SAVEd when needed to be OpenMP shared)
263  INTEGER :: nq
264  INTEGER                       :: nqo, nbtr, nbtr_inca
265  CHARACTER(len=256)              :: type_trac
266  INTEGER,ALLOCATABLE           :: conv_flg(:) ! conv_flg(it)=0 : convection desactivated for tracer number it
267  INTEGER,ALLOCATABLE           :: pbl_flg(:)  ! pbl_flg(it)=0  : boundary layer diffusion desactivaded for tracer number it
268  CHARACTER(len=8),ALLOCATABLE  :: solsym(:)  ! tracer name from inca
269
270  TYPE(t_field),POINTER,SAVE    :: f_ind_cell_glo(:)
271 
272  INTEGER :: iflag_phys   
273
274  INTEGER, ALLOCATABLE, DIMENSION(:) :: hadv_inca  ! index of horizontal trasport schema
275  INTEGER, ALLOCATABLE, DIMENSION(:) :: vadv_inca  ! index of vertical trasport schema
276
277   TYPE(trac_type) ::  tracers_ico2lmdz(nqtot)    !=== TRACERS DESCRIPTORS VECTOR
278   TYPE(isot_type) :: isotopes_ico2lmdz(1)        !=== ISOTOPES PARAMETERS VECTOR
279   INTEGER :: iq
280
281   CHARACTER(LEN=3)      :: descrq(30)            !--- Advection scheme description tags
282   logical, save :: first = .TRUE.
283
284
285    CALL init_distrib_icosa_lmdz
286   
287    ALLOCATE(latfi(klon_omp))
288    ALLOCATE(lonfi(klon_omp))
289    ALLOCATE(airefi(klon_omp))
290    ALLOCATE(bounds_latfi(klon_omp,6))
291    ALLOCATE(bounds_lonfi(klon_omp,6))
292    ALLOCATE(ind_cell_glo_r(klon_omp))
293    ALLOCATE(ind_cell_glo(klon_omp))
294
295    CALL transfer_icosa_to_lmdz(geom%lat_i,latfi)
296    CALL transfer_icosa_to_lmdz(geom%lon_i,lonfi)
297    CALL transfer_icosa_to_lmdz(f_bounds_lat,bounds_latfi)
298    CALL transfer_icosa_to_lmdz(f_bounds_lon,bounds_lonfi)
299    CALL transfer_icosa_to_lmdz(geom%Ai,airefi)
300
301    CALL allocate_field(f_ind_cell_glo,field_t,type_real)
302   
303    DO ind=1,ndomain
304      IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master ) CYCLE
305      CALL swap_dimensions(ind)
306      CALL swap_geometry(ind)
307      DO j=jj_begin,jj_end
308        DO i=ii_begin,ii_end
309          ij=(j-1)*iim+i
310          f_ind_cell_glo(ind)%rval2d(ij)=domain(ind)%assign_cell_glo(i,j)
311        ENDDO
312      ENDDO
313    ENDDO
314
315     
316   CALL transfer_icosa_to_lmdz(f_ind_cell_glo,ind_cell_glo_r)
317   CALL deallocate_field(f_ind_cell_glo)
318   ind_cell_glo=INT(ind_cell_glo_r)
319   DEALLOCATE(ind_cell_glo_r)   
320   
321   
322   CALL reduce_min_lmdz(MINVAL(-ind_cell_glo),ncell_glo_tot) ! reduce_max does not exist in lmdz, use reduce_min
323   CALL bcast_lmdz(ncell_glo_tot)
324   ncell_glo_tot=-ncell_glo_tot
325   ALLOCATE(cell_glo_tot(0:ncell_glo_tot))
326   ALLOCATE(ind_cell_glo_tot(klon_glo))
327   CALL gather_lmdz(ind_cell_glo,ind_cell_glo_tot)
328   CALL bcast_lmdz(ind_cell_glo_tot)
329   
330   cell_glo_tot=-1
331   DO i=1,klon_glo
332     cell_glo_tot(ind_cell_glo_tot(i))= 0
333   ENDDO
334   
335   pos=0
336   DO i=0,ncell_glo_tot
337     IF (cell_glo_tot(i)/=-1) THEN
338       cell_glo_tot(i) = pos
339       pos=pos + 1
340     ENDIF
341   ENDDO
342
343   DO i=1,klon_omp
344     ind_cell_glo(i)=cell_glo_tot(ind_cell_glo(i))
345   ENDDO 
346
347   ind_cell_glo = ind_cell_glo + 1 ! lmdz expect global indices begining to 1 not 0
348 
349   
350!   CALL MPI_ALLGATHER(ncell,1,MPI_INTEGER,ncell_mpi,1,MPI_INTEGER,comm_icosa,ierr)
351!
352!   displ(0)=0
353!   DO i=1,mpi_size-1
354!     displ(i)=displ(i-1)+ncell_mpi(i-1)
355!   ENDDO
356!   
357!   ALLOCATE(ind_glo_tot(ncell_tot))
358!   ALLOCATE(cell_glo_tot(0:ncell_glo_tot-1))
359!   
360!   cell_glo_tot(:)= -1
361!   CALL MPI_ALLGATHERV(ind_glo, ncell, MPI_INTEGER, ind_glo_tot, ncell_mpi, displ, MPI_INTEGER, comm_icosa,ierr)
362!   
363!   DO i=1,ncell_tot
364!     cell_glo_tot(ind_glo_tot(i))= 0
365!   ENDDO
366!
367!   ncell_glo=0
368!   DO i=0,ncell_glo_tot-1
369!     IF (cell_glo_tot(i)/=-1) THEN
370!       cell_glo_tot(i) = ncell_glo
371!       ncell_glo=ncell_glo + 1
372!     ENDIF
373!   ENDDO
374!
375!   DO i=1,ncell
376!     ind_glo(i)=cell_glo_tot(ind_glo(i))
377!   ENDDO 
378
379             
380    CALL init_geometry(klon_omp,lonfi, latfi, bounds_lonfi, bounds_latfi, airefi, ind_cell_glo)
381
382    scaleheight=scale_height/1000. ! Atmospheric scale height (km)
383    aps(1:llm)=0.5*(ap(1:llm)+ap(2:llm+1))
384    bps(1:llm)=0.5*(bp(1:llm)+bp(2:llm+1))
385    pseudoalt(:)=-scaleheight*log(presnivs(:)/preff)
386    CALL init_vertical_layers(llm,preff,scaleheight,ap,bp,aps,bps,presnivs,presinter,pseudoalt)
387
388    ! Initialize tracer names, numbers, etc. for physics
389    !Config  Key  = type_trac
390    !Config  Desc = Choix de couplage avec model de chimie INCA ou REPROBUS
391    !Config  Def  = lmdz
392    !Config  Help =
393    !Config  'lmdz' = pas de couplage, pur LMDZ
394    !Config  'lmdz|inca' = model de chime INCA
395    !Config  'lmdz|repr' = model de chime REPROBUS
396    type_trac = 'lmdz'
397    CALL getin('type_trac',type_trac)
398
399    descrq( 1: 2) = ['LMV','BAK']
400    descrq(10:20) = ['VL1','VLP','FH1','FH2','VLH','   ','PPM','PPS','PPP','   ','SLP']
401    descrq(30)    =  'PRA'
402
403    nqo = 0
404    DO iq=1,nqtot
405
406       tracers_ico2lmdz(iq)%name = tracs(iq)%name
407       
408       tracers_ico2lmdz(iq)%gen0Name = tracs(iq)%name
409       tracers_ico2lmdz(iq)%phase = tracs(iq)%phase
410       tracers_ico2lmdz(iq)%iadv = tracs(iq)%iadv
411         
412       IF (tracs(iq)%component .eq. "dynamico") then
413          tracers_ico2lmdz(iq)%component='lmdz'
414       ELSE
415          tracers_ico2lmdz(iq)%component=tracs(iq)%component
416       ENDIF
417
418       if (tracers_ico2lmdz(iq)%iadv .ne. 0 ) tracers_ico2lmdz(iq)%isAdvected=.true.
419       
420       tracers_ico2lmdz(iq)%longName   = tracers_ico2lmdz(iq)%name
421       IF(tracers_ico2lmdz(iq)%iadv > 0) tracers_ico2lmdz(iq)%longName=TRIM(tracers_ico2lmdz(iq)%name)//descrq(tracers_ico2lmdz(iq)%iadv)
422
423       tracers_ico2lmdz(iq)%isInPhysics= delPhase(tracers_ico2lmdz(iq)%gen0Name) /= 'H2O' .OR. tracers_ico2lmdz(iq)%component /= 'lmdz'
424       tracers_ico2lmdz(iq)%iGeneration = 0
425
426    ENDDO
427
428    nqo = COUNT(delPhase(tracs(:)%name) == 'H2O' .AND. tracers_ico2lmdz(:)%component == 'lmdz') !--- Number of water phases
429
430    isotopes_ico2lmdz(1)%parent='H2O'
431    isotopes_ico2lmdz(1)%phase='gls'
432    isotopes_ico2lmdz(1)%nphas=3
433     
434
435    nbtr=nqtot-nqo
436
437    ALLOCATE(conv_flg(nbtr))
438    ALLOCATE(pbl_flg(nbtr))
439
440    conv_flg(:) = 1 ! convection activated for all tracers
441    pbl_flg(:)  = 1 ! boundary layer activated for all tracers
442       
443
444    CALL init_infotrac_phy
445
446   ! Initialize physical constant
447    day_length=86400
448    CALL getin('day_length',day_length)
449    CALL inifis(day_length,radius,g,kappa*cpp,cpp)
450 
451
452   
453  ! init time
454    annee_ref=2015
455    CALL getin("anneeref",annee_ref)
456   
457    day_ref=1
458    CALL getin("dayref",day_ref)
459   
460    physics_timestep=dt*itau_physics
461    run_length=itaumax*dt
462    ndays=NINT(run_length/day_length)
463   
464    day_ini=INT(itau0*dt/day_length)+day_ref
465    start_time= itau0*dt/day_length-INT(itau0*dt/day_length)
466
467    CALL init_time_lmdz(annee_ref, day_ref, day_ini, start_time, int(ndays), physics_timestep)
468
469!  Additional initializations for aquaplanets
470!    CALL getin("iflag_phys",iflag_phys)
471!    IF (iflag_phys>=100) THEN
472!      CALL iniaqua(klon_omp, iflag_phys)
473!    END IF
474
475 
476
477#ifdef INCA
478    CONTAINS
479
480      SUBROUTINE init_chem_trac()
481        IMPLICIT NONE
482
483        CALL  Init_chem_inca_trac(nbtr)
484
485      END SUBROUTINE init_chem_trac
486
487      SUBROUTINE init_chem_transport()
488
489        IMPLICIT NONE
490
491        CALL init_transport(solsym, conv_flg,pbl_flg, hadv_inca, vadv_inca)
492
493      END SUBROUTINE init_chem_transport
494
495
496#else
497    CONTAINS
498      SUBROUTINE init_chem_trac()
499        IMPLICIT NONE
500
501      END SUBROUTINE init_chem_trac
502
503      SUBROUTINE init_chem_transport()
504
505        IMPLICIT NONE
506
507      END SUBROUTINE init_chem_transport
508
509
510#endif
511
512
513
514  END SUBROUTINE  initialize_physics_omp
515 
516 
517
518
519  SUBROUTINE physics
520  USE icosa
521  USE time_mod
522  USE disvert_mod
523  USE transfert_mod
524  USE mpipara
525  USE xios_mod
526  USE wxios
527  USE trace
528  USE distrib_icosa_lmdz_mod, ONLY : transfer_icosa_to_lmdz, transfer_lmdz_to_icosa
529  USE physics_external_mod, ONLY : it, f_phis, f_ps, f_theta_rhodz, f_u, f_wflux, f_q
530  USE write_field_mod
531  USE checksum_mod
532  USE vorticity_mod
533
534! from LMDZ
535  USE mod_phys_lmdz_omp_data, ONLY: klon_omp
536  USE geometry_mod, ONLY : cell_area
537  USE physiq_mod, ONLY: physiq
538  USE icolmdz_param_gravity_wave, ONLY: param_gravity_wave
539  IMPLICIT NONE
540 
541    REAL(rstd),POINTER :: phis(:)
542    REAL(rstd),POINTER :: ps(:)
543    REAL(rstd),POINTER :: theta_rhodz(:,:,:)
544    REAL(rstd),POINTER :: u(:,:)
545    REAL(rstd),POINTER :: wflux(:,:)
546    REAL(rstd),POINTER :: q(:,:,:)
547    REAL(rstd),POINTER :: p(:,:)
548    REAL(rstd),POINTER :: pks(:)
549    REAL(rstd),POINTER :: pk(:,:)
550    REAL(rstd),POINTER :: p_layer(:,:)
551    REAL(rstd),POINTER :: theta(:,:)
552    REAL(rstd),POINTER :: phi(:,:)
553    REAL(rstd),POINTER :: Temp(:,:)
554    REAL(rstd),POINTER :: ulon(:,:)
555    REAL(rstd),POINTER :: ulat(:,:)
556    REAL(rstd),POINTER :: vort(:,:)
557    REAL(rstd),POINTER :: vortc(:,:)
558    REAL(rstd),POINTER :: dulon(:,:)
559    REAL(rstd),POINTER :: dulat(:,:)
560    REAL(rstd),POINTER :: dTemp(:,:)
561    REAL(rstd),POINTER :: dq(:,:,:)
562    REAL(rstd),POINTER :: dps(:)
563    REAL(rstd),POINTER :: duc(:,:,:)
564
565
566    INTEGER :: ind,l
567   
568    REAL(rstd),ALLOCATABLE,SAVE :: ps_phy(:)
569!$OMP THREADPRIVATE(ps_phy)
570    REAL(rstd),ALLOCATABLE,SAVE :: p_phy(:,:)
571!$OMP THREADPRIVATE(p_phy)
572    REAL(rstd),ALLOCATABLE,SAVE :: p_layer_phy(:,:)
573!$OMP THREADPRIVATE(p_layer_phy)
574    REAL(rstd),ALLOCATABLE,SAVE :: Temp_phy(:,:)
575!$OMP THREADPRIVATE(Temp_phy)
576    REAL(rstd),ALLOCATABLE,SAVE :: phis_phy(:)
577!$OMP THREADPRIVATE(phis_phy)
578    REAL(rstd),ALLOCATABLE,SAVE :: phi_phy(:,:)
579!$OMP THREADPRIVATE(phi_phy)
580    REAL(rstd),ALLOCATABLE,SAVE :: ulon_phy(:,:)
581!$OMP THREADPRIVATE(ulon_phy)
582    REAL(rstd),ALLOCATABLE,SAVE :: ulat_phy(:,:)
583!$OMP THREADPRIVATE(ulat_phy)
584    REAL(rstd),ALLOCATABLE,SAVE :: rot_phy(:,:)
585!$OMP THREADPRIVATE(rot_phy)
586    REAL(rstd),ALLOCATABLE,SAVE :: q_phy(:,:,:)
587!$OMP THREADPRIVATE(q_phy)
588    REAL(rstd),ALLOCATABLE,SAVE :: wflux_phy(:,:)
589!$OMP THREADPRIVATE(wflux_phy)
590    REAL(rstd),ALLOCATABLE,SAVE :: dulon_phy(:,:)
591!$OMP THREADPRIVATE(dulon_phy)
592    REAL(rstd),ALLOCATABLE,SAVE :: dulat_phy(:,:)
593!$OMP THREADPRIVATE(dulat_phy)
594    REAL(rstd),ALLOCATABLE,SAVE :: dTemp_phy(:,:)
595!$OMP THREADPRIVATE(dTemp_phy)
596    REAL(rstd),ALLOCATABLE,SAVE :: dq_phy(:,:,:)
597!$OMP THREADPRIVATE(dq_phy)
598    REAL(rstd),ALLOCATABLE,SAVE :: dps_phy(:)
599!$OMP THREADPRIVATE(dps_phy)
600    REAL(rstd)   :: dtphy
601    LOGICAL      :: debut
602    LOGICAL      :: lafin
603    LOGICAL,SAVE :: first=.TRUE.
604!$OMP THREADPRIVATE(first)
605
606   
607    IF(first) THEN
608      debut=.TRUE.
609    ELSE
610      debut=.FALSE.
611    ENDIF
612
613
614    IF(it-itau0>=itaumax) THEN
615      lafin=.TRUE.
616    ELSE
617      lafin=.FALSE.
618    ENDIF
619
620    IF (first) THEN
621      first=.FALSE.
622      CALL init_message(f_u,req_e1_vect,req_u)
623      CALL init_message(f_vort,req_z1_scal,req_z)
624      ALLOCATE(ps_phy(klon_omp))
625      ALLOCATE(p_phy(klon_omp,llm+1))
626      ALLOCATE(p_layer_phy(klon_omp,llm))
627      ALLOCATE(Temp_phy(klon_omp,llm))
628      ALLOCATE(phis_phy(klon_omp))
629      ALLOCATE(phi_phy(klon_omp,llm))
630      ALLOCATE(ulon_phy(klon_omp,llm))
631      ALLOCATE(ulat_phy(klon_omp,llm))
632      ALLOCATE(rot_phy(klon_omp,llm))
633      ALLOCATE(q_phy(klon_omp,llm,nqtot))
634      ALLOCATE(wflux_phy(klon_omp,llm))
635      ALLOCATE(dulon_phy(klon_omp,llm))
636      ALLOCATE(dulat_phy(klon_omp,llm))
637      ALLOCATE(dTemp_phy(klon_omp,llm))
638      ALLOCATE(dq_phy(klon_omp,llm,nqtot))
639      ALLOCATE(dps_phy(klon_omp))
640!$OMP BARRIER
641     
642!      CALL param_gravity_wave
643    ENDIF
644
645
646!$OMP MASTER       
647!    CALL update_calendar(it)
648!$OMP END MASTER
649!$OMP BARRIER
650    dtphy=itau_physics*dt
651   
652   
653   
654    CALL transfert_message(f_u,req_u)
655    DO ind=1,ndomain
656      IF (assigned_domain(ind)) THEN
657        CALL swap_dimensions(ind)
658        CALL swap_geometry(ind)
659        u=f_u(ind)
660        vort=f_vort(ind)
661        CALL compute_vorticity(u,vort)
662      ENDIF
663    ENDDO
664
665    CALL transfert_message(f_vort,req_z)
666
667   
668    DO ind=1,ndomain
669      CALL swap_dimensions(ind)
670      IF (assigned_domain(ind)) THEN
671        CALL swap_geometry(ind)
672     
673        phis=f_phis(ind)
674        ps=f_ps(ind)
675        theta_rhodz=f_theta_rhodz(ind)
676        u=f_u(ind)
677        q=f_q(ind)
678        wflux=f_wflux(ind)
679        p=f_p(ind)
680        pks=f_pks(ind)
681        pk=f_pk(ind)
682        p_layer=f_p_layer(ind)
683        theta=f_theta(ind)
684        phi=f_phi(ind)
685        Temp=f_Temp(ind)
686        ulon=f_ulon(ind)
687        ulat=f_ulat(ind)
688        vort=f_vort(ind)
689        vortc=f_vortc(ind)
690           
691        CALL grid_icosa_to_physics
692
693      ENDIF
694    ENDDO
695   
696!$OMP BARRIER
697!$OMP MASTER
698    CALL SYSTEM_CLOCK(start_clock)
699!$OMP END MASTER
700    CALL trace_start("physic")
701!    CALL trace_off()
702
703
704!    CALL writeField("p_in",f_p)
705!    CALL writeField("p_layer_in",f_p_layer)
706!    CALL writeField("phi_in",f_phi)
707!    CALL writeField("phis_in",f_phis)
708!    CALL writeField("ulon_in",f_ulon)
709!    CALL writeField("ulat_in",f_ulat)
710!    CALL writeField("Temp_in",f_Temp)
711!    CALL writeField("q_in",f_q)
712!    CALL writeField("wflux_in",f_wflux)
713!     CALL writeField("vortc",f_vortc)
714
715!    CALL checksum(f_p)
716!    CALL checksum(f_p_layer)
717!    CALL checksum(f_phi)
718!    CALL checksum(f_phis)
719!    CALL checksum(f_ulon)
720!    CALL checksum(f_ulat)
721!    CALL checksum(f_Temp)
722!    CALL checksum(f_q)
723!    CALL checksum(f_wflux)
724
725    CALL transfer_icosa_to_lmdz(f_p      , p_phy)
726    CALL transfer_icosa_to_lmdz(f_p_layer, p_layer_phy)
727    CALL transfer_icosa_to_lmdz(f_phi    , phi_phy)
728    CALL transfer_icosa_to_lmdz(f_phis   , phis_phy )
729    CALL transfer_icosa_to_lmdz(f_ulon   , ulon_phy )
730    CALL transfer_icosa_to_lmdz(f_ulat   , ulat_phy)
731    CALL transfer_icosa_to_lmdz(f_vortc   , rot_phy)
732    CALL transfer_icosa_to_lmdz(f_Temp   , Temp_phy)
733    CALL transfer_icosa_to_lmdz(f_q      , q_phy)
734    CALL transfer_icosa_to_lmdz(f_wflux  , wflux_phy)
735
736    DO l=1,llm
737      wflux_phy(:,l) = - wflux_phy(:,l)*cell_area(:)
738      phi_phy(:,l)=phi_phy(:,l)-phis_phy(:)
739    ENDDO
740   
741    CALL wxios_set_context()
742 
743    ! Ehouarn: rot_phy() not implemented!! Set it to zero for now
744!    rot_phy(:,:)=0
745    CALL physiq(klon_omp, llm, debut, lafin, dtphy, &
746                p_phy, p_layer_phy, phi_phy, phis_phy, presnivs, &
747                ulon_phy, ulat_phy, rot_phy, Temp_phy, q_phy, wflux_phy, &
748                dulon_phy, dulat_phy, dTemp_phy, dq_phy, dps_phy)
749   
750    CALL transfer_lmdz_to_icosa(dulon_phy, f_dulon )
751    CALL transfer_lmdz_to_icosa(dulat_phy, f_dulat )
752    CALL transfer_lmdz_to_icosa(dTemp_phy, f_dTemp )
753    CALL transfer_lmdz_to_icosa(dq_phy   , f_dq )
754    CALL transfer_lmdz_to_icosa(dps_phy  , f_dps )
755 
756!    CALL writeField("dulon_out",f_dulon)
757!    CALL writeField("dulat_out",f_dulat)
758!    CALL writeField("dTemp_out",f_dTemp)
759!    CALL writeField("dq_out",f_dq)
760!    CALL writeField("dps_out",f_dps)
761
762!    CALL checksum(f_dulon)
763!    CALL checksum(f_dulat)
764!    CALL checksum(f_dTemp)
765!    CALL checksum(f_dq)
766!    CALL checksum(f_dps)
767   
768    CALL send_message(f_dps,req_dps0)
769    CALL send_message(f_dulon,req_dulon0)
770    CALL send_message(f_dulat,req_dulat0)
771    CALL send_message(f_dTemp,req_dTemp0)
772    CALL send_message(f_dq,req_dq0)
773
774    CALL wait_message(req_dps0)
775    CALL wait_message(req_dulon0)
776    CALL wait_message(req_dulat0)
777    CALL wait_message(req_dTemp0)
778    CALL wait_message(req_dq0)
779
780
781!    CALL trace_on()
782    CALL trace_end("physic")
783!$OMP MASTER
784    CALL SYSTEM_CLOCK(stop_clock)
785    count_clock=count_clock+stop_clock-start_clock
786!$OMP END MASTER
787
788!$OMP BARRIER                       
789
790    DO ind=1,ndomain
791      CALL swap_dimensions(ind)
792      IF (assigned_domain(ind)) THEN
793        CALL swap_geometry(ind)
794
795        theta_rhodz=f_theta_rhodz(ind)
796        u=f_u(ind)
797        q=f_q(ind)
798        ps=f_ps(ind)
799        dulon=f_dulon(ind)
800        dulat=f_dulat(ind)
801        Temp=f_temp(ind)
802        dTemp=f_dTemp(ind)
803        dq=f_dq(ind)
804        dps=f_dps(ind)
805        duc=f_duc(ind)
806        p=f_p(ind)
807        pks=f_pks(ind)
808        pk=f_pk(ind)
809     
810        CALL grid_physics_to_icosa
811      ENDIF
812    ENDDO
813
814!$OMP BARRIER
815    CALL xios_set_context   
816   
817 
818  CONTAINS
819
820    SUBROUTINE grid_icosa_to_physics
821    USE pression_mod
822    USE exner_mod
823    USE theta2theta_rhodz_mod
824    USE geopotential_mod
825    USE wind_from_lonlat_mod
826    USE omp_para
827    IMPLICIT NONE
828   
829    REAL(rstd) :: uc(3)
830    INTEGER :: i,j,ij,l
831
832! compute pression
833
834      DO    l    = ll_begin,ll_endp1
835        DO j=jj_begin,jj_end
836          DO i=ii_begin,ii_end
837            ij=(j-1)*iim+i
838            p(ij,l) = ap(l) + bp(l) * ps(ij)
839          ENDDO
840        ENDDO
841      ENDDO
842
843!$OMP BARRIER
844
845! compute exner
846       
847       IF (is_omp_first_level) THEN
848         DO j=jj_begin,jj_end
849            DO i=ii_begin,ii_end
850               ij=(j-1)*iim+i
851               pks(ij) = cpp * ( ps(ij)/preff ) ** kappa
852            ENDDO
853         ENDDO
854       ENDIF
855
856       ! 3D : pk
857       DO l = ll_begin,ll_end
858          DO j=jj_begin,jj_end
859             DO i=ii_begin,ii_end
860                ij=(j-1)*iim+i
861                pk(ij,l) = cpp * ((.5/preff)*(p(ij,l)+p(ij,l+1))) ** kappa
862             ENDDO
863          ENDDO
864       ENDDO
865
866!$OMP BARRIER
867
868!   compute theta, temperature and pression at layer
869    DO    l    = ll_begin, ll_end
870      DO j=jj_begin,jj_end
871        DO i=ii_begin,ii_end
872          ij=(j-1)*iim+i
873          theta(ij,l) = theta_rhodz(ij,l,1) / ((p(ij,l)-p(ij,l+1))/g)
874          Temp(ij,l) = theta(ij,l) * pk(ij,l) / cpp
875          p_layer(ij,l)=preff*(pk(ij,l)/cpp)**(1./kappa)
876        ENDDO
877      ENDDO
878    ENDDO
879
880
881!!! Compute geopotential
882       
883  ! for first layer
884  IF (is_omp_first_level) THEN
885    DO j=jj_begin,jj_end
886      DO i=ii_begin,ii_end
887        ij=(j-1)*iim+i
888        phi( ij,1 ) = phis( ij ) + theta(ij,1) * ( pks(ij) - pk(ij,1) )
889      ENDDO
890    ENDDO
891  ENDIF
892!!-> implicit flush on phi(:,1)
893
894!$OMP BARRIER
895         
896  ! for other layers
897  DO l = ll_beginp1, ll_end
898    DO j=jj_begin,jj_end
899      DO i=ii_begin,ii_end
900        ij=(j-1)*iim+i
901        phi(ij,l) =  0.5 * ( theta(ij,l)  + theta(ij,l-1) )  &
902                         * (  pk(ij,l-1) -  pk(ij,l)    )
903      ENDDO
904    ENDDO
905  ENDDO       
906
907!$OMP BARRIER
908
909
910  IF (is_omp_first_level) THEN
911    DO l = 2, llm
912      DO j=jj_begin,jj_end
913! ---> Bug compilo intel ici en openmp
914! ---> Couper la boucle
915       IF (j==jj_end+1) PRINT*,"this message must not be printed"
916        DO i=ii_begin,ii_end
917          ij=(j-1)*iim+i
918          phi(ij,l) = phi(ij,l)+ phi(ij,l-1)
919        ENDDO
920      ENDDO
921    ENDDO
922! --> IMPLICIT FLUSH on phi --> non
923  ENDIF
924
925! compute wind centered lon lat compound
926    DO l=ll_begin,ll_end
927      DO j=jj_begin,jj_end
928        DO i=ii_begin,ii_end
929          ij=(j-1)*iim+i
930          uc(:)=1/Ai(ij)*                                                                                                &
931                        ( 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,:))  &
932                         + 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,:))          &
933                         + 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,:))          &
934                         + 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,:))    &
935                         + 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,:))&
936                         + 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,:)))
937          ulon(ij,l)=sum(uc(:)*elon_i(ij,:))
938          ulat(ij,l)=sum(uc(:)*elat_i(ij,:))
939        ENDDO
940      ENDDO
941    ENDDO
942
943
944! compute centered vorticity
945   
946    DO l=ll_begin,ll_end
947      DO j=jj_begin,jj_end
948        DO i=ii_begin,ii_end
949          ij=(j-1)*iim+i
950          vortc(ij,l) =  Riv(ij,vup)    * vort(ij+z_up,l)    + &
951                         Riv(ij,vlup)  * vort(ij+z_lup,l)   + &
952                         Riv(ij,vldown)* vort(ij+z_ldown,l) + &
953                         Riv(ij,vdown) * vort(ij+z_down,l)  + &
954                         Riv(ij,vrdown)* vort(ij+z_rdown,l) + &
955                         Riv(ij,vrup)  * vort(ij+z_rup,l)
956      ENDDO
957    ENDDO
958  ENDDO
959
960
961
962!$OMP BARRIER
963    END SUBROUTINE grid_icosa_to_physics
964
965
966    SUBROUTINE grid_physics_to_icosa
967    USE theta2theta_rhodz_mod
968    USE omp_para
969    IMPLICIT NONE
970      INTEGER :: i,j,ij,l,iq
971         
972      DO l=ll_begin,ll_end
973        DO j=jj_begin,jj_end
974          DO i=ii_begin,ii_end
975            ij=(j-1)*iim+i
976            duc(ij,:,l)=dulon(ij,l)*elon_i(ij,:)+dulat(ij,l)*elat_i(ij,:)
977          ENDDO
978        ENDDO
979      ENDDO
980
981      DO l=ll_begin,ll_end
982        DO j=jj_begin,jj_end
983          DO i=ii_begin,ii_end
984            ij=(j-1)*iim+i
985            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,:) )
986            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,:) )
987            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,:) )
988          ENDDO
989        ENDDO
990      ENDDO         
991
992      DO l=ll_begin,ll_end
993        DO j=jj_begin,jj_end
994          DO i=ii_begin,ii_end
995            ij=(j-1)*iim+i
996            Temp(ij,l)=Temp(ij,l)+ dtphy * dTemp(ij,l)
997          ENDDO
998        ENDDO
999      ENDDO         
1000     
1001      DO iq=1,nqtot
1002        DO l=ll_begin,ll_end
1003          DO j=jj_begin,jj_end
1004            DO i=ii_begin,ii_end
1005              ij=(j-1)*iim+i
1006              q(ij,l,iq)=q(ij,l,iq)+ dtphy * dq(ij,l,iq)
1007            ENDDO
1008          ENDDO
1009        ENDDO
1010      ENDDO
1011
1012!$OMP BARRIER
1013     
1014     IF (is_omp_first_level) THEN
1015       DO j=jj_begin,jj_end
1016         DO i=ii_begin,ii_end
1017           ij=(j-1)*iim+i
1018           ps(ij)=ps(ij)+ dtphy * dps(ij)
1019          ENDDO
1020       ENDDO
1021     ENDIF
1022
1023!     CALL compute_temperature2theta_rhodz(ps,Temp,theta_rhodz,0)
1024
1025! compute pression
1026!$OMP BARRIER
1027      DO    l    = ll_begin,ll_endp1
1028        DO j=jj_begin,jj_end
1029          DO i=ii_begin,ii_end
1030            ij=(j-1)*iim+i
1031            p(ij,l) = ap(l) + bp(l) * ps(ij)
1032          ENDDO
1033        ENDDO
1034      ENDDO
1035
1036!$OMP BARRIER
1037
1038! compute exner
1039       
1040       IF (is_omp_first_level) THEN
1041         DO j=jj_begin,jj_end
1042            DO i=ii_begin,ii_end
1043               ij=(j-1)*iim+i
1044               pks(ij) = cpp * ( ps(ij)/preff ) ** kappa
1045            ENDDO
1046         ENDDO
1047       ENDIF
1048
1049       ! 3D : pk
1050       DO l = ll_begin,ll_end
1051          DO j=jj_begin,jj_end
1052             DO i=ii_begin,ii_end
1053                ij=(j-1)*iim+i
1054                pk(ij,l) = cpp * ((.5/preff)*(p(ij,l)+p(ij,l+1))) ** kappa
1055             ENDDO
1056          ENDDO
1057       ENDDO
1058
1059!$OMP BARRIER
1060
1061!   compute theta, temperature and pression at layer
1062    DO    l    = ll_begin, ll_end
1063      DO j=jj_begin,jj_end
1064        DO i=ii_begin,ii_end
1065          ij=(j-1)*iim+i
1066          theta_rhodz(ij,l,1) = temp(ij,l) * ((p(ij,l)-p(ij,l+1))/g) / (pk(ij,l) / cpp )
1067        ENDDO
1068      ENDDO
1069    ENDDO
1070   
1071    END SUBROUTINE grid_physics_to_icosa
1072
1073
1074
1075  END SUBROUTINE physics
1076
1077
1078
1079
1080
1081END MODULE interface_icosa_lmdz_mod
Note: See TracBrowser for help on using the repository browser.