source: LMDZ5/trunk/libf/dyn3dmem/call_calfis_mod.F90 @ 2375

Last change on this file since 2375 was 2375, checked in by Ehouarn Millour, 9 years ago

Fix in the computation of the date; the convention is that it corresponds to the time at the end of the current dynamics or physics step (except when in backward Matsuno step where it remains unchanged as it was correctly updated during the forward part of the step).
Note that the relationship between itau and date is a bit tricky as itau is incremented between Matsuno forward and backward steps (and from a leapfrog step to the next) but unchanged from Matsuno bacward step to leapfrog step.
Because of change in date when calling physics, bench results will change with this revision.
EM

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