LMDZ
call_calfis_mod.F90
Go to the documentation of this file.
1 !#define DEBUG_IO
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 
32 CONTAINS
33 
34  SUBROUTINE call_calfis_allocate
35  USE bands
37  USE parallel_lmdz
38  USE dimensions_mod
39  USE infotrac
40  IMPLICIT NONE
41  TYPE(distrib),POINTER :: d
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)
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
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
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, &
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
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)
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
343  else
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)
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)
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
409  !$OMP END MASTER
410 
411  END SUBROUTINE call_calfis
412 
413 END MODULE call_calfis_mod
Definition: bands.F90:4
subroutine top_bound_loc(vcov, ucov, teta, masse, dt)
Definition: top_bound_loc.F:5
real, dimension(:), pointer, save dpfi
!$Header llmm1 INTEGER ip1jmp1
Definition: paramet.h:14
real, dimension(:,:,:), pointer, save q
real, dimension(:,:,:), pointer, save dq
subroutine exner_milieu_loc(ngrid, ps, p, pks, pk, pkf)
!$Header llmp1
Definition: paramet.h:14
!$Id bp(llm+1)
Definition: vampir.F90:1
real, dimension(:,:), pointer, save dufi
integer, save ij_end
subroutine vtb(number)
Definition: vampir.F90:52
subroutine massdair_loc(p, masse)
Definition: massdair_loc.F:2
subroutine calfis_loc(lafin, jD_cur, jH_cur, pucov, pvcov, pteta, pq, pmasse, pps, pp, ppk, pphis, pphi, pducov, pdvcov, pdteta, pdq, flxw, pdufi, pdvfi, pdhfi, pdqfi, pdpsfi)
Definition: calfis_loc.F:28
character(len=10), save planet_type
Definition: control_mod.F90:32
real, dimension(:,:,:), pointer, save dqfi
!$Id Turb_fcg_gcssold get_uvd hqturb_gcssold endif!large scale llm day day1 day day1 *dt_toga endif!time annee_ref dt_toga u_toga vq_toga w_prof vq_prof llm day day1 day day1 *dt_dice endif!time annee_ref dt_dice swup_dice vg_dice omega_dice tg_prof vg_profd w_profd omega_profd!do llm!print llm l llm
integer, save day_step
Definition: control_mod.F90:15
real, dimension(:,:), pointer, save masse
integer, parameter timer_caldyn
Definition: times.F90:7
subroutine resume_timer(no_timer)
Definition: times.F90:87
subroutine exner_hyb_loc(ngrid, ps, p, pks, pk, pkf)
integer, save nqtot
Definition: infotrac.F90:6
real, dimension(:,:), pointer, save dvfi
!$Id ysinus ok_strato
Definition: logic.h:10
real, dimension(:,:), pointer, save teta
integer, save ijb_v
subroutine pression_loc(ngrid, ap, bp, ps, p)
Definition: pression_loc.F:2
!$Id && day_ini
Definition: temps.h:15
!$Id mode_top_bound COMMON comconstr dtphys
Definition: comconst.h:7
real, dimension(:), pointer, save phis
type(distrib), target, save distrib_physic
Definition: bands.F90:21
logical, save pole_nord
!$Id day_ref
Definition: temps.h:15
real, dimension(:,:), pointer, save phi
real, dimension(:,:), pointer, save vcov
real, dimension(:,:), pointer, save dteta
!FH On elimine toutes les clefs physiques dans la dynamique prt_level
integer, parameter vthallo
Definition: vampir.F90:7
subroutine register_hallo_u(Field, ll, RUp, Rdown, SUp, SDown, a_request)
Definition: mod_hallo.F90:942
Definition: times.F90:1
subroutine set_distrib(d)
subroutine sendrequest(a_Request)
Definition: mod_hallo.F90:1072
real, dimension(:,:), pointer, save pkf
subroutine call_calfis_allocate
integer, save ij_begin
integer, save ije_v
real, dimension(:,:), pointer, save flxw
subroutine call_calfis(itau, lafin, ucov_dyn, vcov_dyn, teta_dyn, masse_dyn, ps_dyn, phis_dyn, q_dyn, flxw_dyn)
subroutine vte(number)
Definition: vampir.F90:69
real, dimension(:,:), pointer, save pk
real, dimension(:,:), pointer, save dv
real, dimension(:,:), pointer, save du
!$Id forward
Definition: logic.h:10
subroutine suspend_timer(no_timer)
Definition: times.F90:70
!$Id leapf
Definition: logic.h:10
c c zjulian c cym CALL iim cym klev iim
Definition: ini_bilKP_ave.h:24
integer, save ip_ebil_dyn
Definition: control_mod.F90:29
real, dimension(:), pointer, save pks
character(len=maxlen) function int2str(int)
subroutine geopot_loc(ngrid, teta, pk, pks, phis, phi)
Definition: geopot_loc.F:2
subroutine settag(a_request, tag)
Definition: mod_hallo.F90:180
type(distrib), target, save distrib_physic_bis
Definition: bands.F90:22
integer, save ije_u
real, dimension(:,:), pointer, save p
real, dimension(:), pointer, save ps
real, dimension(:,:), pointer, save ucov
type(distrib), target, save distrib_caldyn
Definition: bands.F90:17
!$Id start_time
Definition: temps.h:15
subroutine register_hallo_v(Field, ll, RUp, Rdown, SUp, SDown, a_request)
Definition: mod_hallo.F90:1007
subroutine waitrequest(a_Request)
Definition: mod_hallo.F90:1196
real, dimension(:,:), pointer, save dtetafi
integer, save ijnb_u
subroutine addfi_loc(pdt, leapf, forward, pucov, pvcov, pteta, pq, pps, pdufi, pdvfi, pdhfi, pdqfi, pdpfi)
Definition: addfi_loc.F:7
!$Header!gestion des impressions de sorties et de débogage la sortie standard prt_level COMMON comprint lunout
Definition: iniprint.h:7
integer, save ijb_u
integer, parameter vtphysiq
Definition: vampir.F90:8