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

Last change on this file since 1985 was 1982, checked in by emillour, 7 years ago

Dynamico interface:
Follow-up of having the possibility of tracer names up to 30 characters long.
EM

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