source: trunk/ICOSA_LMDZ/src/phymars/interface_icosa_lmdz.f90 @ 2425

Last change on this file since 2425 was 2424, checked in by adelavois, 5 years ago

Update of DYNAMICO / Mars interface with isotopes nqperes and nqfils

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