My Project
 All Classes Files Functions Variables Macros
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 :: alpha(:,:)
15  REAL,POINTER,SAVE :: beta(:,:)
16  REAL,POINTER,SAVE :: pks(:)
17  REAL,POINTER,SAVE :: pk(:,:)
18  REAL,POINTER,SAVE :: pkf(:,:)
19  REAL,POINTER,SAVE :: phi(:,:)
20  REAL,POINTER,SAVE :: du(:,:)
21  REAL,POINTER,SAVE :: dv(:,:)
22  REAL,POINTER,SAVE :: dteta(:,:)
23  REAL,POINTER,SAVE :: dq(:,:,:)
24  REAL,POINTER,SAVE :: dufi(:,:)
25  REAL,POINTER,SAVE :: dvfi(:,:)
26  REAL,POINTER,SAVE :: dtetafi(:,:)
27  REAL,POINTER,SAVE :: dqfi(:,:,:)
28  REAL,POINTER,SAVE :: dpfi(:)
29 
30 
31 
32 
33 
34 CONTAINS
35 
37  USE bands
38  USE allocate_field
39  USE parallel
40  USE dimensions
41  USE infotrac
42  IMPLICIT NONE
43  TYPE(distrib),POINTER :: d
44  d=>distrib_physic
45 
46  CALL allocate_u(ucov,llm,d)
47  CALL allocate_v(vcov,llm,d)
48  CALL allocate_u(teta,llm,d)
49  CALL allocate_u(masse,llm,d)
50  CALL allocate_u(ps,d)
51  CALL allocate_u(phis,d)
52  CALL allocate_u(q,llm,nqtot,d)
53  CALL allocate_u(flxw,llm,d)
54  CALL allocate_u(p,llmp1,d)
55  CALL allocate_u(alpha,llm,d)
56  CALL allocate_u(beta,llm,d)
57  CALL allocate_u(pks,d)
58  CALL allocate_u(pk,llm,d)
59  CALL allocate_u(pkf,llm,d)
60  CALL allocate_u(phi,llm,d)
61  CALL allocate_u(du,llm,d)
62  CALL allocate_v(dv,llm,d)
63  CALL allocate_u(dteta,llm,d)
64  CALL allocate_u(dq,llm,nqtot,d)
65  CALL allocate_u(dufi,llm,d)
66  CALL allocate_v(dvfi,llm,d)
67  CALL allocate_u(dtetafi,llm,d)
68  CALL allocate_u(dqfi,llm,nqtot,d)
69  CALL allocate_u(dpfi,d)
70 
71  END SUBROUTINE call_calfis_allocate
72 
73 
74  SUBROUTINE call_calfis(itau,lafin,clesphy0,ucov_dyn,vcov_dyn,teta_dyn,masse_dyn,ps_dyn, &
75  phis_dyn,q_dyn,flxw_dyn)
76  USE dimensions
77  USE parallel
78  USE times
79  USE mod_hallo
80  USE bands
81  USE vampir
82  USE infotrac
83  USE control_mod
84  USE write_field_loc
85  USE write_field
86  IMPLICIT NONE
87  include "comconst.h"
88  include "comvert.h"
89  include "logic.h"
90  include "temps.h"
91  include "iniprint.h"
92 
93  REAL :: clesphy0( : )
94  INTEGER :: itau
95  LOGICAL :: lafin
96  REAL :: ucov_dyn(ijb_u:ije_u,llm)
97  REAL :: vcov_dyn(ijb_v:ije_v,llm)
98  REAL :: teta_dyn(ijb_u:ije_u,llm)
99  REAL :: masse_dyn(ijb_u:ije_u,llm)
100  REAL :: ps_dyn(ijb_u:ije_u)
101  REAL :: phis_dyn(ijb_u:ije_u)
102  REAL :: q_dyn(ijb_u:ije_u,llm,nqtot)
103  REAL :: flxw_dyn(ijb_u:ije_u,llm)
104 
105  REAL :: dufi_tmp(iip1,llm)
106  REAL :: dvfi_tmp(iip1,llm)
107  REAL :: dtetafi_tmp(iip1,llm)
108  REAL :: dpfi_tmp(iip1)
109  REAL :: dqfi_tmp(iip1,llm,nqtot)
110 
111  REAL :: jd_cur, jh_cur
112  CHARACTER(LEN=15) :: ztit
113  TYPE(request) :: request_physic
114  INTEGER :: ijb,ije,l,j
115 
116 
117 #ifdef DEBUG_IO
118  CALL writefield_u('ucovfi',ucov)
119  CALL writefield_v('vcovfi',vcov)
120  CALL writefield_u('tetafi',teta)
121  CALL writefield_u('pfi',p)
122  CALL writefield_u('pkfi',pk)
123  DO j=1,nqtot
124  CALL writefield_u('qfi'//trim(int2str(j)),q(:,:,j))
125  ENDDO
126 #endif
127 
128 !
129 ! ....... Ajout P.Le Van ( 17/04/96 ) ...........
130 !
131 
132 
133  !$OMP MASTER
134  CALL suspend_timer(timer_caldyn)
135  WRITE(lunout,*) 'leapfrog_p: Entree dans la physique : Iteration No ',itau
136  !$OMP END MASTER
137 
138  jd_cur = jd_ref + day_ini - day_ref &
139  & + itau/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,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,alpha,beta, 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  CALL calfis_loc(lafin ,jd_cur, jh_cur, &
230  ucov,vcov,teta,q,masse,ps,p,pk,phis,phi , &
231  du,dv,dteta,dq, &
232  flxw, &
233  clesphy0, dufi,dvfi,dtetafi,dqfi,dpfi )
234 
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  IF (ok_strato) THEN
324  CALL top_bound_loc( vcov,ucov,teta,masse,dufi,dvfi,dtetafi)
325  ENDIF
326 
327 #ifdef DEBUG_IO
328  CALL writefield_u('ucovfi',ucov)
329  CALL writefield_v('vcovfi',vcov)
330  CALL writefield_u('tetafi',teta)
331  CALL writefield_u('psfi',ps)
332  DO j=1,nqtot
333  CALL writefield_u('qfi'//trim(int2str(j)),q(:,:,j))
334  ENDDO
335 #endif
336 
337  CALL addfi_loc( dtphys, leapf, forward , &
338  ucov, vcov, teta , q ,ps , &
339  dufi, dvfi, dtetafi , dqfi ,dpfi )
340 
341 #ifdef DEBUG_IO
342  CALL writefield_u('ucovfi',ucov)
343  CALL writefield_v('vcovfi',vcov)
344  CALL writefield_u('tetafi',teta)
345  CALL writefield_u('psfi',ps)
346  DO j=1,nqtot
347  CALL writefield_u('qfi'//trim(int2str(j)),q(:,:,j))
348  ENDDO
349 #endif
350 
351  !$OMP BARRIER
352  !$OMP MASTER
353  CALL vte(vtphysiq)
354  CALL vtb(vthallo)
355  !$OMP END MASTER
356 
357  CALL settag(request_physic,800)
358  CALL register_swapfield_u(ucov,ucov_dyn,distrib_caldyn,request_physic)
359  CALL register_swapfield_v(vcov,vcov_dyn,distrib_caldyn,request_physic)
360  CALL register_swapfield_u(teta,teta_dyn,distrib_caldyn,request_physic)
361  CALL register_swapfield_u(masse,masse_dyn,distrib_caldyn,request_physic)
362  CALL register_swapfield_u(ps,ps_dyn,distrib_caldyn,request_physic)
363  CALL register_swapfield_u(q,q_dyn,distrib_caldyn,request_physic)
364  CALL sendrequest(request_physic)
365  !$OMP BARRIER
366  CALL waitrequest(request_physic)
367 
368  !$OMP BARRIER
369  !$OMP MASTER
370  CALL vte(vthallo)
371  CALL set_distrib(distrib_caldyn)
372  !$OMP END MASTER
373  !$OMP BARRIER
374 
375 !
376 ! Diagnostique de conservation de l'energie : difference
377  IF (ip_ebil_dyn.ge.1 ) THEN
378  ztit='bil phys'
379  CALL diagedyn(ztit,2,1,1,dtphys,ucov, vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
380  ENDIF
381 
382 #ifdef DEBUG_IO
383  CALL writefield_u('ucovfi',ucov_dyn)
384  CALL writefield_v('vcovfi',vcov_dyn)
385  CALL writefield_u('tetafi',teta_dyn)
386  CALL writefield_u('psfi',ps_dyn)
387  DO j=1,nqtot
388  CALL writefield_u('qfi'//trim(int2str(j)),q_dyn(:,:,j))
389  ENDDO
390 #endif
391 
392 
393 !-jld
394  !$OMP MASTER
395  CALL resume_timer(timer_caldyn)
396  !$OMP END MASTER
397 
398  END SUBROUTINE call_calfis
399 
400 END MODULE call_calfis_mod