source: LMDZ6/trunk/libf/dyn3dmem/lmdz_call_calfis.F90 @ 5072

Last change on this file since 5072 was 5066, checked in by abarral, 5 months ago

Transform gr_dyn_fi_p.F, gr_fi_dyn_p.F, calfis_loc.F into free-form modules.
Reorder CPP_PARA keys in lmdz_call_calfis.F90, lmdz_calfis_loc.F90, lmdz_gr_dyn_fi_p.F90, lmdz_gr_fi_dyn_p.F90 to avoid implicit declarations.
Remove redundant -cpp -D.. on arch.
Correct "!OMP" -> "!$OMP"
Correct typo in lmdz_xios.F90, wstats.F90

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 12.1 KB
Line 
1!#define DEBUG_IO
2MODULE lmdz_call_calfis
3
4    REAL,POINTER,SAVE :: ucov(:,:)
5    REAL,POINTER,SAVE :: vcov(:,:)
6    REAL,POINTER,SAVE :: teta(:,:)
7    REAL,POINTER,SAVE :: masse(:,:)
8    REAL,POINTER,SAVE :: ps(:)
9    REAL,POINTER,SAVE :: phis(:)
10    REAL,POINTER,SAVE :: q(:,:,:)
11    REAL,POINTER,SAVE :: flxw(:,:)
12
13    REAL,POINTER,SAVE :: p(:,:)
14    REAL,POINTER,SAVE :: pks(:)
15    REAL,POINTER,SAVE :: pk(:,:)
16    REAL,POINTER,SAVE :: pkf(:,:)
17    REAL,POINTER,SAVE :: phi(:,:)
18    REAL,POINTER,SAVE :: du(:,:)
19    REAL,POINTER,SAVE :: dv(:,:)
20    REAL,POINTER,SAVE :: dteta(:,:)
21    REAL,POINTER,SAVE :: dq(:,:,:)
22    REAL,POINTER,SAVE :: dufi(:,:)
23    REAL,POINTER,SAVE :: dvfi(:,:)
24    REAL,POINTER,SAVE :: dtetafi(:,:)
25    REAL,POINTER,SAVE :: dqfi(:,:,:)
26    REAL,POINTER,SAVE :: dpfi(:)
27   
28   
29   
30   
31   
32CONTAINS
33
34  SUBROUTINE call_calfis_allocate
35  USE bands
36  USE allocate_field_mod
37  USE parallel_lmdz
38  USE dimensions_mod
39  USE infotrac, ONLY: nqtot
40  IMPLICIT NONE
41    TYPE(distrib),POINTER :: d
42    d=>distrib_physic
43
44    CALL allocate_u(ucov,llm,d)
45    CALL allocate_v(vcov,llm,d)
46    CALL allocate_u(teta,llm,d)
47    CALL allocate_u(masse,llm,d)
48    CALL allocate_u(ps,d)
49    CALL allocate_u(phis,d)
50    CALL allocate_u(q,llm,nqtot,d)
51    CALL allocate_u(flxw,llm,d)
52    CALL allocate_u(p,llmp1,d)
53    CALL allocate_u(pks,d)
54    pks(:)=0
55    CALL allocate_u(pk,llm,d)
56    pk(:,:)=0
57    CALL allocate_u(pkf,llm,d)
58    CALL allocate_u(phi,llm,d)
59    CALL allocate_u(du,llm,d)
60    CALL allocate_v(dv,llm,d)
61    CALL allocate_u(dteta,llm,d)
62    CALL allocate_u(dq,llm,nqtot,d)
63    CALL allocate_u(dufi,llm,d)
64    CALL allocate_v(dvfi,llm,d)
65    CALL allocate_u(dtetafi,llm,d)
66    CALL allocate_u(dqfi,llm,nqtot,d)
67    CALL allocate_u(dpfi,d)
68 
69  END SUBROUTINE call_calfis_allocate
70 
71 
72  SUBROUTINE call_calfis(itau,lafin,ucov_dyn,vcov_dyn,teta_dyn,masse_dyn,ps_dyn, &
73                         phis_dyn,q_dyn,flxw_dyn)
74  USE dimensions_mod
75  use exner_hyb_loc_m, only: exner_hyb_loc
76  use exner_milieu_loc_m, only: exner_milieu_loc
77  USE parallel_lmdz
78  USE times
79  USE mod_hallo
80  USE Bands
81  USE vampir
82  USE infotrac, ONLY: nqtot
83  USE control_mod
84  USE write_field_loc
85  USE write_field
86  USE comconst_mod, ONLY: dtphys
87  USE logic_mod, ONLY: leapf, forward, ok_strato
88  USE comvert_mod, ONLY: ap, bp, pressure_exner
89  USE temps_mod, ONLY: day_ini, day_ref, jd_ref, jh_ref, start_time
90#ifdef CPP_PHYS
91  USE lmdz_calfis_loc
92#endif
93
94  IMPLICIT NONE
95    INCLUDE "iniprint.h"
96
97    INTEGER,INTENT(IN) :: itau ! (time) iteration step number
98    LOGICAL,INTENT(IN) :: lafin ! .true. if final time step
99    REAL,INTENT(INOUT) :: ucov_dyn(ijb_u:ije_u,llm) ! covariant zonal wind
100    REAL,INTENT(INOUT) :: vcov_dyn(ijb_v:ije_v,llm) ! covariant meridional wind
101    REAL,INTENT(INOUT) :: teta_dyn(ijb_u:ije_u,llm) ! potential temperature
102    REAL,INTENT(INOUT) :: masse_dyn(ijb_u:ije_u,llm) ! air mass
103    REAL,INTENT(INOUT) :: ps_dyn(ijb_u:ije_u) ! surface pressure
104    REAL,INTENT(INOUT) :: phis_dyn(ijb_u:ije_u) ! surface geopotential
105    REAL,INTENT(INOUT) :: q_dyn(ijb_u:ije_u,llm,nqtot) ! advected tracers
106    REAL,INTENT(INOUT) :: flxw_dyn(ijb_u:ije_u,llm) ! vertical mass flux
107
108    REAL :: dufi_tmp(iip1,llm)   
109    REAL :: dvfi_tmp(iip1,llm) 
110    REAL :: dtetafi_tmp(iip1,llm)
111    REAL :: dpfi_tmp(iip1)
112    REAL :: dqfi_tmp(iip1,llm,nqtot)
113
114    REAL :: jD_cur, jH_cur
115    CHARACTER(LEN=15) :: ztit
116    TYPE(Request),SAVE :: Request_physic
117!$OMP THREADPRIVATE(Request_physic )
118    INTEGER :: ijb,ije,l,iq
119   
120   
121#ifdef DEBUG_IO   
122    CALL WriteField_u('ucovfi',ucov)
123    CALL WriteField_v('vcovfi',vcov)
124    CALL WriteField_u('tetafi',teta)
125    CALL WriteField_u('pfi',p)
126    CALL WriteField_u('pkfi',pk)
127    DO iq=1,nqtot
128      CALL WriteField_u('qfi'//trim(int2str(iq)),q(:,:,iq))
129    ENDDO
130#endif
131
132!
133!     .......   Ajout   P.Le Van ( 17/04/96 )   ...........
134!
135
136
137  !$OMP MASTER
138    CALL suspend_timer(timer_caldyn)
139    IF (prt_level >= 10) THEN
140      WRITE(lunout,*) 'leapfrog_p: Entree dans la physique : Iteration No ',itau
141    ENDIF
142  !$OMP END MASTER
143   
144           jD_cur = jD_ref + day_ini - day_ref                           &
145     &        + (itau+1)/day_step
146
147           IF (planet_type .eq."generic") THEN
148              ! AS: we make jD_cur to be pday
149              jD_cur = int(day_ini + itau/day_step)
150           ENDIF
151
152           jH_cur = jH_ref + start_time +                                &
153     &              mod(itau+1,day_step)/float(day_step)
154    if (jH_cur > 1.0 ) then
155      jD_cur = jD_cur +1.
156      jH_cur = jH_cur -1.
157    endif
158
159!   Inbterface avec les routines de phylmd (phymars ... )
160!   -----------------------------------------------------
161
162!+jld
163
164!  Diagnostique de conservation de l'energie : initialisation
165 
166!-jld
167  !$OMP BARRIER
168  !$OMP MASTER
169    CALL VTb(VThallo)
170  !$OMP END MASTER
171
172#ifdef DEBUG_IO   
173    CALL WriteField_u('ucovfi',ucov)
174    CALL WriteField_v('vcovfi',vcov)
175    CALL WriteField_u('tetafi',teta)
176    CALL WriteField_u('pfi',p)
177    CALL WriteField_u('pkfi',pk)
178#endif
179   
180    CALL SetTag(Request_physic,800)
181    CALL Register_SwapField_u(ucov_dyn,ucov,distrib_physic,Request_physic,up=2,down=2)
182    CALL Register_SwapField_v(vcov_dyn,vcov,distrib_physic,Request_physic,up=2,down=2)
183    CALL Register_SwapField_u(teta_dyn,teta,distrib_physic,Request_physic,up=2,down=2)
184    CALL Register_SwapField_u(masse_dyn,masse,distrib_physic,Request_physic,up=1,down=2)
185    CALL Register_SwapField_u(ps_dyn,ps,distrib_physic,Request_physic,up=2,down=2)
186    CALL Register_SwapField_u(phis_dyn,phis,distrib_physic,Request_physic,up=2,down=2)
187    CALL Register_SwapField_u(q_dyn,q,distrib_physic,Request_physic,up=2,down=2)
188    CALL Register_SwapField_u(flxw_dyn,flxw,distrib_physic,Request_physic,up=2,down=2)
189 
190    CALL SendRequest(Request_Physic)
191  !$OMP BARRIER
192    CALL WaitRequest(Request_Physic)       
193
194  !$OMP BARRIER
195  !$OMP MASTER
196    CALL Set_Distrib(distrib_Physic)
197    CALL VTe(VThallo)
198       
199    CALL VTb(VTphysiq)
200  !$OMP END MASTER
201  !$OMP BARRIER
202
203    CALL pression_loc (  ip1jmp1, ap, bp, ps,  p      )
204
205  !$OMP BARRIER
206    CALL exner_hyb_loc(  ip1jmp1, ps, p, pks, pk, pkf )
207  !$OMP BARRIER
208    CALL geopot_loc  ( ip1jmp1, teta  , pk , pks,  phis  , phi   )
209
210
211    CALL Register_Hallo_u(p,llmp1,2,2,2,2,Request_physic)
212    CALL Register_Hallo_u(pk,llm,2,2,2,2,Request_physic)
213    CALL Register_Hallo_u(phi,llm,2,2,2,2,Request_physic)
214       
215    CALL SendRequest(Request_Physic)
216  !$OMP BARRIER
217    CALL WaitRequest(Request_Physic)
218             
219  !$OMP BARRIER
220 
221 
222#ifdef DEBUG_IO   
223    CALL WriteField_u('ucovfi',ucov)
224    CALL WriteField_v('vcovfi',vcov)
225    CALL WriteField_u('tetafi',teta)
226    CALL WriteField_u('pfi',p)
227    CALL WriteField_u('pkfi',pk)
228    DO iq=1,nqtot
229      CALL WriteField_u('qfi'//trim(int2str(iq)),q(:,:,iq))
230    ENDDO
231#endif
232
233  !$OMP BARRIER
234
235#ifdef CPP_PHYS
236    CALL calfis_loc(lafin ,jD_cur, jH_cur,                       &
237                     ucov,vcov,teta,q,masse,ps,p,pk,phis,phi ,   &
238                     du,dv,dteta,dq,                             &
239                     flxw, dufi,dvfi,dtetafi,dqfi,dpfi  )
240#endif
241    ijb=ij_begin
242    ije=ij_end 
243    IF ( .not. pole_nord) THEN
244 
245    !$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
246      DO l=1,llm
247        dufi_tmp(1:iip1,l)   = dufi(ijb:ijb+iim,l)
248        dvfi_tmp(1:iip1,l)   = dvfi(ijb:ijb+iim,l) 
249        dtetafi_tmp(1:iip1,l)= dtetafi(ijb:ijb+iim,l) 
250        dqfi_tmp(1:iip1,l,:) = dqfi(ijb:ijb+iim,l,:) 
251      ENDDO
252    !$OMP END DO NOWAIT
253
254    !$OMP MASTER
255      dpfi_tmp(1:iip1)     = dpfi(ijb:ijb+iim) 
256    !$OMP END MASTER
257   
258    ENDIF ! of if ( .not. pole_nord)
259
260  !$OMP BARRIER
261  !$OMP MASTER
262    CALL Set_Distrib(distrib_Physic_bis)
263    CALL VTb(VThallo)
264  !$OMP END MASTER
265  !$OMP BARRIER
266 
267    CALL Register_Hallo_u(dufi,llm,1,0,0,1,Request_physic)
268    CALL Register_Hallo_v(dvfi,llm,1,0,0,1,Request_physic)
269    CALL Register_Hallo_u(dtetafi,llm,1,0,0,1,Request_physic)
270    CALL Register_Hallo_u(dpfi,1,1,0,0,1,Request_physic)
271
272    DO iq=1,nqtot
273      CALL Register_Hallo_u(dqfi(:,:,iq),llm,1,0,0,1,Request_physic)
274    ENDDO
275       
276    CALL SendRequest(Request_Physic)
277  !$OMP BARRIER
278    CALL WaitRequest(Request_Physic)
279             
280  !$OMP BARRIER
281  !$OMP MASTER
282    CALL VTe(VThallo)
283    CALL Set_Distrib(distrib_Physic)
284  !$OMP END MASTER
285  !$OMP BARRIER       
286    ijb=ij_begin
287    IF (.not. pole_nord) THEN
288       
289    !$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
290      DO l=1,llm
291        dufi(ijb:ijb+iim,l) = dufi(ijb:ijb+iim,l)+dufi_tmp(1:iip1,l)
292        dvfi(ijb:ijb+iim,l) = dvfi(ijb:ijb+iim,l)+dvfi_tmp(1:iip1,l)
293        dtetafi(ijb:ijb+iim,l) = dtetafi(ijb:ijb+iim,l)+dtetafi_tmp(1:iip1,l)
294        dqfi(ijb:ijb+iim,l,:) = dqfi(ijb:ijb+iim,l,:) + dqfi_tmp(1:iip1,l,:)
295      ENDDO
296    !$OMP END DO NOWAIT
297
298    !$OMP MASTER
299      dpfi(ijb:ijb+iim)   = dpfi(ijb:ijb+iim)+ dpfi_tmp(1:iip1)
300    !$OMP END MASTER
301         
302    endif ! of if (.not. pole_nord)
303       
304       
305#ifdef DEBUG_IO           
306    CALL WriteField_u('dufi',dufi)
307    CALL WriteField_v('dvfi',dvfi)
308    CALL WriteField_u('dtetafi',dtetafi)
309    CALL WriteField_u('dpfi',dpfi)
310    DO iq=1,nqtot
311      CALL WriteField_u('dqfi'//trim(int2str(iq)),dqfi(:,:,iq))
312    ENDDO
313#endif
314
315  !$OMP BARRIER
316
317!      ajout des tendances physiques:
318!      ------------------------------
319#ifdef DEBUG_IO   
320    CALL WriteField_u('ucovfi',ucov)
321    CALL WriteField_v('vcovfi',vcov)
322    CALL WriteField_u('tetafi',teta)
323    CALL WriteField_u('psfi',ps)
324    DO iq=1,nqtot
325      CALL WriteField_u('qfi'//trim(int2str(iq)),q(:,:,iq))
326    ENDDO
327#endif
328
329#ifdef DEBUG_IO           
330    CALL WriteField_u('ucovfi',ucov)
331    CALL WriteField_v('vcovfi',vcov)
332    CALL WriteField_u('tetafi',teta)
333    CALL WriteField_u('psfi',ps)
334    DO iq=1,nqtot
335      CALL WriteField_u('qfi'//trim(int2str(iq)),q(:,:,iq))
336    ENDDO
337#endif
338
339    CALL addfi_loc( dtphys, leapf, forward   ,              &
340                    ucov, vcov, teta , q   ,ps ,            &
341                    dufi, dvfi, dtetafi , dqfi ,dpfi  )
342    ! since addfi updates ps(), also update p(), masse() and pk()
343    CALL pression_loc(ip1jmp1,ap,bp,ps,p)
344!$OMP BARRIER
345    CALL massdair_loc(p,masse)
346!$OMP BARRIER
347    if (pressure_exner) then
348      CALL exner_hyb_loc(ijnb_u,ps,p,pks,pk,pkf)
349    else
350      CALL exner_milieu_loc(ijnb_u,ps,p,pks,pk,pkf)
351    endif
352!$OMP BARRIER
353
354#ifdef DEBUG_IO   
355    CALL WriteField_u('ucovfi',ucov)
356    CALL WriteField_v('vcovfi',vcov)
357    CALL WriteField_u('tetafi',teta)
358    CALL WriteField_u('psfi',ps)
359    DO iq=1,nqtot
360      CALL WriteField_u('qfi'//trim(int2str(iq)),q(:,:,iq))
361    ENDDO
362#endif
363
364    IF (ok_strato) THEN
365!      CALL top_bound_loc( vcov,ucov,teta,masse,dufi,dvfi,dtetafi)
366      CALL top_bound_loc(vcov,ucov,teta,masse,dtphys)
367    ENDIF
368
369  !$OMP BARRIER
370  !$OMP MASTER
371    CALL VTe(VTphysiq)
372    CALL VTb(VThallo)
373  !$OMP END MASTER
374
375    CALL SetTag(Request_physic,800)
376    CALL Register_SwapField_u(ucov,ucov_dyn,distrib_caldyn,Request_physic)
377    CALL Register_SwapField_v(vcov,vcov_dyn,distrib_caldyn,Request_physic)
378    CALL Register_SwapField_u(teta,teta_dyn,distrib_caldyn,Request_physic)
379    CALL Register_SwapField_u(masse,masse_dyn,distrib_caldyn,Request_physic)
380    CALL Register_SwapField_u(ps,ps_dyn,distrib_caldyn,Request_physic)
381    CALL Register_SwapField_u(q,q_dyn,distrib_caldyn,Request_physic)
382    CALL SendRequest(Request_Physic)
383  !$OMP BARRIER
384    CALL WaitRequest(Request_Physic)     
385
386  !$OMP BARRIER
387  !$OMP MASTER
388    CALL VTe(VThallo)
389    CALL set_distrib(distrib_caldyn)
390  !$OMP END MASTER
391  !$OMP BARRIER
392
393!
394!  Diagnostique de conservation de l'energie : difference
395    IF (ip_ebil_dyn.ge.1 ) THEN
396      ztit='bil phys'
397!      CALL diagedyn(ztit,2,1,1,dtphys,ucov, vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
398      write(lunout,*)"call_calfis: diagedyn disabled in dyn3dmem !!"
399    ENDIF
400
401#ifdef DEBUG_IO   
402    CALL WriteField_u('ucovfi',ucov_dyn)
403    CALL WriteField_v('vcovfi',vcov_dyn)
404    CALL WriteField_u('tetafi',teta_dyn)
405    CALL WriteField_u('psfi',ps_dyn)
406    DO iq=1,nqtot
407      CALL WriteField_u('qfi'//trim(int2str(iq)),q_dyn(:,:,iq))
408    ENDDO
409#endif
410
411
412!-jld
413    !$OMP MASTER
414      CALL resume_timer(timer_caldyn)
415    !$OMP END MASTER
416
417  END SUBROUTINE call_calfis
418 
419END MODULE lmdz_call_calfis
Note: See TracBrowser for help on using the repository browser.