source: trunk/ICOSA_LMDZ/src/phyvenus/interface_icosa_lmdz.f90 @ 2368

Last change on this file since 2368 was 2330, checked in by slebonnois, 5 years ago

SL: modifications to be able to turn the physics with Cp(T)

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