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

Last change on this file since 3655 was 2603, checked in by Ehouarn Millour, 8 years ago

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