source: LMDZ5/trunk/libf/dyn3dmem/call_dissip_mod.F90 @ 1848

Last change on this file since 1848 was 1848, checked in by yann meurdesoif, 11 years ago

Solve performance problem comming from declarations of a derived type.

YM

File size: 8.6 KB
Line 
1MODULE call_dissip_mod
2
3    REAL,POINTER,SAVE :: ucov(:,:)
4    REAL,POINTER,SAVE :: vcov(:,:)
5    REAL,POINTER,SAVE :: teta(:,:)
6    REAL,POINTER,SAVE :: p(:,: )
7    REAL,POINTER,SAVE :: pk(:,:)
8
9    REAL,POINTER,SAVE :: ucont(:,:)
10    REAL,POINTER,SAVE :: vcont(:,:)
11    REAL,POINTER,SAVE :: ecin(:,:)
12    REAL,POINTER,SAVE :: ecin0(:,:)
13    REAL,POINTER,SAVE :: dudis(:,:)
14    REAL,POINTER,SAVE :: dvdis(:,:)
15    REAL,POINTER,SAVE :: dtetadis(:,:)
16    REAL,POINTER,SAVE :: dtetaecdt(:,:)
17
18
19
20CONTAINS
21 
22  SUBROUTINE call_dissip_allocate
23  USE bands
24  USE allocate_field_mod
25  USE parallel_lmdz
26  USE dimensions_mod
27  USE dissip_mod, ONLY : dissip_allocate
28  IMPLICIT NONE
29    TYPE(distrib),POINTER :: d
30    d=>distrib_dissip
31
32    CALL allocate_u(ucov,llm,d)
33    CALL allocate_v(vcov,llm,d)
34    CALL allocate_u(teta,llm,d)
35    CALL allocate_u(p,llmp1,d)
36    CALL allocate_u(pk,llm,d)
37    CALL allocate_u(ucont,llm,d)
38    CALL allocate_v(vcont,llm,d)
39    CALL allocate_u(ecin,llm,d)
40    CALL allocate_u(ecin0,llm,d)
41    CALL allocate_u(dudis,llm,d)
42    CALL allocate_v(dvdis,llm,d)
43    CALL allocate_u(dtetadis,llm,d)
44    CALL allocate_u(dtetaecdt,llm,d)
45   
46   
47    CALL dissip_allocate
48   
49  END SUBROUTINE call_dissip_allocate
50 
51  SUBROUTINE call_dissip_switch_dissip(dist)
52  USE allocate_field_mod
53  USE bands
54  USE parallel_lmdz
55  USE dissip_mod, ONLY : dissip_switch_dissip
56  IMPLICIT NONE
57    TYPE(distrib),INTENT(IN) :: dist
58
59    CALL switch_u(ucov,distrib_dissip,dist)
60    CALL switch_v(vcov,distrib_dissip,dist)
61    CALL switch_u(teta,distrib_dissip,dist)
62    CALL switch_u(p,distrib_dissip,dist)
63    CALL switch_u(pk,distrib_dissip,dist)
64    CALL switch_u(ucont,distrib_dissip,dist)
65    CALL switch_v(vcont,distrib_dissip,dist)
66    CALL switch_u(ecin,distrib_dissip,dist)
67    CALL switch_u(ecin0,distrib_dissip,dist)
68    CALL switch_u(dudis,distrib_dissip,dist)
69    CALL switch_v(dvdis,distrib_dissip,dist)
70    CALL switch_u(dtetadis,distrib_dissip,dist)
71    CALL switch_u(dtetaecdt,distrib_dissip,dist)
72
73    CALL dissip_switch_dissip(dist)
74   
75  END SUBROUTINE call_dissip_switch_dissip 
76 
77
78 
79  SUBROUTINE call_dissip(ucov_dyn,vcov_dyn,teta_dyn,p_dyn,pk_dyn,ps_dyn)
80  USE dimensions_mod
81  USE parallel_lmdz
82  USE times
83  USE mod_hallo
84  USE Bands
85  USE vampir
86  USE write_field_loc
87  IMPLICIT NONE
88    INCLUDE 'comgeom.h'
89    REAL :: ucov_dyn(ijb_u:ije_u,llm)
90    REAL :: vcov_dyn(ijb_v:ije_v,llm)
91    REAL :: teta_dyn(ijb_u:ije_u,llm)
92    REAL :: p_dyn(ijb_u:ije_u,llmp1 )
93    REAL :: pk_dyn(ijb_u:ije_u,llm)
94    REAL :: ps_dyn(ijb_u:ije_u)
95    REAL :: tppn(iim),tpps(iim)
96    REAL :: tpn,tps
97
98    REAL  SSUM
99    LOGICAL,PARAMETER :: dissip_conservative=.TRUE.
100    TYPE(Request),SAVE :: Request_dissip
101!$OMP THREADPRIVATE(Request_dissip )   
102    INTEGER :: ij,l,ijb,ije
103 
104   
105  !$OMP MASTER
106    CALL suspend_timer(timer_caldyn)
107       
108!       print*,'Entree dans la dissipation : Iteration No ',true_itau
109!   calcul de l'energie cinetique avant dissipation
110!       print *,'Passage dans la dissipation'
111
112    CALL VTb(VThallo)
113  !$OMP END MASTER
114
115  !$OMP BARRIER
116
117    CALL Register_SwapField_u(ucov_dyn,ucov,distrib_dissip, Request_dissip,up=1,down=1)
118    CALL Register_SwapField_v(vcov_dyn,vcov,distrib_dissip, Request_dissip,up=1,down=1)
119    CALL Register_SwapField_u(teta_dyn,teta,distrib_dissip, Request_dissip)
120    CALL Register_SwapField_u(p_dyn,p,distrib_dissip,Request_dissip)
121    CALL Register_SwapField_u(pk_dyn,pk,distrib_dissip,Request_dissip)
122
123    CALL SendRequest(Request_dissip)       
124  !$OMP BARRIER
125    CALL WaitRequest(Request_dissip)       
126
127  !$OMP BARRIER
128  !$OMP MASTER
129    CALL set_distrib(distrib_dissip)
130    CALL VTe(VThallo)
131    CALL VTb(VTdissipation)
132    CALL start_timer(timer_dissip)
133  !$OMP END MASTER
134  !$OMP BARRIER
135
136    CALL covcont_loc(llm,ucov,vcov,ucont,vcont)
137    CALL enercin_loc(vcov,ucov,vcont,ucont,ecin0)
138
139!   dissipation
140
141!        CALL FTRACE_REGION_BEGIN("dissip")
142    CALL dissip_loc(vcov,ucov,teta,p,dvdis,dudis,dtetadis)
143
144#ifdef DEBUG_IO   
145    CALL WriteField_u('dudis',dudis)
146    CALL WriteField_v('dvdis',dvdis)
147    CALL WriteField_u('dtetadis',dtetadis)
148#endif
149 
150!      CALL FTRACE_REGION_END("dissip")
151         
152    ijb=ij_begin
153    ije=ij_end
154  !$OMP DO SCHEDULE(STATIC,OMP_CHUNK)       
155    DO l=1,llm
156      ucov(ijb:ije,l)=ucov(ijb:ije,l)+dudis(ijb:ije,l)
157    ENDDO
158  !$OMP END DO NOWAIT       
159
160    IF (pole_sud) ije=ije-iip1
161   
162  !$OMP DO SCHEDULE(STATIC,OMP_CHUNK)       
163    DO l=1,llm
164      vcov(ijb:ije,l)=vcov(ijb:ije,l)+dvdis(ijb:ije,l)
165    ENDDO
166  !$OMP END DO NOWAIT       
167
168!       teta=teta+dtetadis
169
170
171!------------------------------------------------------------------------
172    IF (dissip_conservative) THEN
173!       On rajoute la tendance due a la transform. Ec -> E therm. cree
174!       lors de la dissipation
175    !$OMP BARRIER
176    !$OMP MASTER
177      CALL suspend_timer(timer_dissip)
178      CALL VTb(VThallo)
179    !$OMP END MASTER
180      CALL Register_Hallo_u(ucov,llm,1,1,1,1,Request_Dissip)
181      CALL Register_Hallo_v(vcov,llm,1,1,1,1,Request_Dissip)
182      CALL SendRequest(Request_Dissip)
183    !$OMP BARRIER
184      CALL WaitRequest(Request_Dissip)
185    !$OMP MASTER
186      CALL VTe(VThallo)
187      CALL resume_timer(timer_dissip)
188    !$OMP END MASTER
189    !$OMP BARRIER           
190      CALL covcont_loc(llm,ucov,vcov,ucont,vcont)
191      CALL enercin_loc(vcov,ucov,vcont,ucont,ecin)
192           
193      ijb=ij_begin
194      ije=ij_end
195    !$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
196      DO l=1,llm
197        DO ij=ijb,ije
198           dtetaecdt(ij,l)= (ecin0(ij,l)-ecin(ij,l))/ pk(ij,l)
199           dtetadis(ij,l)=dtetadis(ij,l)+dtetaecdt(ij,l)
200        ENDDO
201      ENDDO
202    !$OMP END DO NOWAIT           
203
204    ENDIF
205
206    ijb=ij_begin
207    ije=ij_end
208
209  !$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
210    DO l=1,llm
211      DO ij=ijb,ije
212         teta(ij,l)=teta(ij,l)+dtetadis(ij,l)
213      ENDDO
214    ENDDO
215  !$OMP END DO NOWAIT         
216
217!------------------------------------------------------------------------
218
219
220!    .......        P. Le Van (  ajout  le 17/04/96  )   ...........
221!   ...      Calcul de la valeur moyenne, unique de h aux poles  .....
222!
223
224    ijb=ij_begin
225    ije=ij_end
226         
227    IF (pole_nord) THEN
228 
229   !$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
230      DO l  =  1, llm
231        DO ij =  1,iim
232          tppn(ij)  = aire(  ij    ) * teta(  ij    ,l)
233        ENDDO
234        tpn  = SSUM(iim,tppn,1)/apoln
235
236        DO ij = 1, iip1
237          teta(  ij    ,l) = tpn
238        ENDDO
239      ENDDO
240    !$OMP END DO NOWAIT
241
242         if (1 == 0) then
243!!! Ehouarn: lines here 1) kill 1+1=2 in the dynamics
244!!!                     2) should probably not be here anyway
245!!! but are kept for those who would want to revert to previous behaviour
246    !$OMP MASTER               
247      DO ij =  1,iim
248        tppn(ij)  = aire(  ij    ) * ps_dyn (  ij    )
249      ENDDO
250      tpn  = SSUM(iim,tppn,1)/apoln
251 
252      DO ij = 1, iip1
253        ps_dyn(  ij    ) = tpn
254      ENDDO
255    !$OMP END MASTER
256   
257    ENDIF ! of if (1 == 0)
258    endif ! of of (pole_nord)
259       
260    IF (pole_sud) THEN
261
262    !$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
263      DO l  =  1, llm
264        DO ij =  1,iim
265          tpps(ij)  = aire(ij+ip1jm) * teta(ij+ip1jm,l)
266        ENDDO
267       
268        tps  = SSUM(iim,tpps,1)/apols
269
270        DO ij = 1, iip1
271          teta(ij+ip1jm,l) = tps
272        ENDDO
273      ENDDO
274    !$OMP END DO NOWAIT
275
276    if (1 == 0) then
277!!! Ehouarn: lines here 1) kill 1+1=2 in the dynamics
278!!!                     2) should probably not be here anyway
279!!! but are kept for those who would want to revert to previous behaviour
280    !$OMP MASTER               
281      DO ij =  1,iim
282        tpps(ij)  = aire(ij+ip1jm) * ps_dyn (ij+ip1jm)
283      ENDDO
284      tps  = SSUM(iim,tpps,1)/apols
285 
286      DO ij = 1, iip1
287        ps_dyn(ij+ip1jm) = tps
288      ENDDO
289    !$OMP END MASTER
290    ENDIF ! of if (1 == 0)
291    endif ! of if (pole_sud)
292
293
294  !$OMP BARRIER
295  !$OMP MASTER
296    CALL VTe(VTdissipation)
297    CALL stop_timer(timer_dissip)
298    CALL VTb(VThallo)
299  !$OMP END MASTER
300 
301    CALL Register_SwapField_u(ucov,ucov_dyn,distrib_caldyn,Request_dissip)
302    CALL Register_SwapField_v(vcov,vcov_dyn,distrib_caldyn,Request_dissip)
303    CALL Register_SwapField_u(teta,teta_dyn,distrib_caldyn,Request_dissip)
304    CALL Register_SwapField_u(p,p_dyn,distrib_caldyn,Request_dissip)
305    CALL Register_SwapField_u(pk,pk_dyn,distrib_caldyn,Request_dissip)
306
307    CALL SendRequest(Request_dissip)       
308
309  !$OMP BARRIER
310    CALL WaitRequest(Request_dissip)       
311  !$OMP BARRIER
312  !$OMP MASTER
313    CALL set_distrib(distrib_caldyn)
314    CALL VTe(VThallo)
315    CALL resume_timer(timer_caldyn)
316!        print *,'fin dissipation'
317  !$OMP END MASTER
318  !$OMP BARRIER
319 
320 
321  END SUBROUTINE call_dissip
322
323END MODULE call_dissip_mod
Note: See TracBrowser for help on using the repository browser.