source: LMDZ6/trunk/libf/dyn3dmem/call_calfis_mod.F90 @ 5217

Last change on this file since 5217 was 5084, checked in by Laurent Fairhead, 4 months ago

Reverting to r4065. Updating fortran standard broke too much stuff. Will do it by smaller chunks
AB, LF

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