source: trunk/SOURCES/resol_adv_diff_2D-sept2009.f90

Last change on this file was 285, checked in by aquiquet, 5 years ago

Small modifications for gfortran enabled compilations

File size: 11.0 KB
Line 
1!> \file resol_adv_diff_2D-sept2009.f90
2!!Resoud l equation adv-diffusion par une methode de relaxation
3!<
4
5!> \namespace reso_adv_diff_2D_vect
6!! Resoud l equation adv-diffusion par une methode de relaxation
7!! @note   - passage du vecteur et de l epaisseur  en dummy
8!! @note   - epaisseur prescrite selon un masque
9!! @note use module3D_phy
10!! @todo essayer dans le futur la methode de Richard dUX/dx = H dU/dx + U dH/dx
11!<
12
13module reso_adv_diff_2D_vect
14
15use module3D_phy
16
17
18
19implicit none
20double precision :: omega   !< parametre schema temporel de la resolution partie diffusion
21double precision :: mu_adv  !< parametre schema temporel de la resolution advection
22double precision :: upwind  !< schema spatial pour l'advection
23
24! tableaux de travail. resolution M H = Frelax
25double precision,dimension(nx,ny) :: crelax      !< diagnonale de M
26double precision,dimension(nx,ny) :: arelax      !< sous diagonale selon x
27double precision,dimension(nx,ny) :: brelax      !< sur diagonale selon x
28double precision,dimension(nx,ny) :: drelax      !< sous diagonale selon y
29double precision,dimension(nx,ny) :: erelax      !< sur diagonale selon y
30double precision,dimension(nx,ny) :: frelax      !< vecteur
31double precision,dimension(nx,ny) :: c_west      !< sur demi mailles Ux
32double precision,dimension(nx,ny) :: c_east      !< sur demi mailles Ux
33double precision,dimension(nx,ny) :: c_north     !< sur demi mailles Uy
34double precision,dimension(nx,ny) :: c_south     !< sur demi mailles Uy
35
36double precision,dimension(nx,ny) :: bdx         !< pente socle
37double precision,dimension(nx,ny) :: bdy         !< pente socle
38
39double precision,dimension(nx,ny) :: hdx         !< pente epaisseur
40double precision,dimension(nx,ny) :: hdy         !< pente epaisseur
41
42
43
44
45integer, dimension(2) :: ijmax    ! position de maxval
46integer :: iFAIL
47
48contains 
49
50!> initialise  init_reso_adv_diff_2D
51!! definition des parametres qui gerent la resolution
52!! @return omega, mu_adv, upwind
53
54subroutine init_reso_adv_diff_2D
55
56
57!write(num_rep_42,*)'partie diffusive'
58!write(num_rep_42,*)'le type de schema temporel diffusif depend de omega'
59!write(num_rep_42,*)'0 -> explicite, 0.5 -> Crank-Nicolson'
60!write(num_rep_42,*)'1 -> semi implicite,  >1 -> over-implicite'
61
62! over implicite propose par Hindmasrsh
63omega = 2.5
64!write(num_rep_42,*)'omega = ',omega
65!write(num_rep_42,*)
66!write(num_rep_42,*)'partie advective'
67!write(num_rep_42,*)' le schéma temporel advectif dépend de mu'
68!write(num_rep_42,*)'0 -> explicite, 0.5 -> Crank-Nicolson'
69!write(num_rep_42,*)'1 -> semi implicite,  >1 -> over-implicite'
70mu_adv=1.  ! l'over implicite ne s'applique pas a l'equation advective
71!write(num_rep_42,*)'mu_adv = ',mu_adv
72!write(num_rep_42,*)' le schéma spatial dépend de upwind'
73!write(num_rep_42,*)'1 -> schema upwind,  0.5 -> schema centre'
74upwind=1
75!write(num_rep_42,*)'upwind = ',upwind
76!write(num_rep_42,*)'------------------------------------------------------'
77
78return
79end subroutine init_reso_adv_diff_2D
80!------------------------------------------------------------------
81
82!> subroutine resol_adv_diff_2D_Hpresc
83!! definition des parametres qui gerent la resolution
84!! @param Dfx,Dfy             termes diffusifs
85!! @param advx,advy           termes advectifs
86!! @param vieuxH              ancienne valeur de H
87!! @param H_presc             la valeur H prescrit pour certains noeuds
88!! @param i_Hpresc            le masque ou la valeur de H est prescrite
89!! @param bil                 le bilan pour la colonne
90!! @return newH               la valeur de H apres le pas de temps
91
92subroutine resol_adv_diff_2D_vect(Dfx,Dfy,advx,advy,H_presc,i_Hpresc,bil,vieuxH,newH)
93
94!$ USE OMP_LIB
95implicit none
96
97double precision,dimension(nx,ny), intent(in) :: Dfx          !< terme diffusif selon x
98double precision,dimension(nx,ny), intent(in) :: Dfy          !< terme diffusif selon y
99double precision,dimension(nx,ny), intent(in) :: Advx         !< terme advectif selon x
100double precision,dimension(nx,ny), intent(in) :: Advy         !< terme advectif selon y
101double precision,dimension(nx,ny), intent(in) :: vieuxH       !< ancienne valeur de H
102double precision,dimension(nx,ny), intent(out):: newH         !< nouvelle valeur de H
103
104double precision,dimension(nx,ny), intent(in) :: bil          !< bilan de masse pour la colonne
105
106double precision,dimension(nx,ny), intent(in) :: H_presc      !< H value if prescribed
107integer,dimension(nx,ny), intent(in) :: i_Hpresc  !< 1 if H is prescribed on this node, else 0
108
109
110! tableaux de travail. resolution M H = Frelax
111!!$real,dimension(nx,ny) :: crelax      !< diagnonale de M
112!!$real,dimension(nx,ny) :: arelax      !< sous diagonale selon x
113!!$real,dimension(nx,ny) :: brelax      !< sur diagonale selon x
114!!$real,dimension(nx,ny) :: drelax      !< sous diagonale selon y
115!!$real,dimension(nx,ny) :: erelax      !< sur diagonale selon y
116!!$real,dimension(nx,ny) :: frelax      !< vecteur
117!!$real,dimension(nx,ny) :: c_west      !< sur demi mailles Ux
118!!$real,dimension(nx,ny) :: c_east      !< sur demi mailles Ux
119!!$real,dimension(nx,ny) :: c_north     !< sur demi mailles Uy
120!!$real,dimension(nx,ny) :: c_south     !<sur demi mailles Ux
121
122!!$real,dimension(nx,ny) :: bdx         !< pente socle
123!!$real,dimension(nx,ny) :: bdy         !< pente socle
124!!$
125!!$real,dimension(nx,ny) :: hdx         !< pente epaisseur
126!!$real,dimension(nx,ny) :: hdy         !< pente epaisseur
127
128double precision :: frdx,frdy                    !< pour calcul frelax : termes diffusion
129double precision :: fraxw,fraxe,frays,frayn      !< termes advection
130
131double precision,dimension(nx,ny) :: deltah      ! dans calcul relax
132double precision :: delh                         ! dans calcul relax
133double precision :: testh                        ! dans calcul relax
134
135logical :: stopp
136integer :: ntour
137double precision :: reste
138
139
140if (itracebug.eq.1)  call tracebug(' Entree dans routine resolution_diffusion')
141
142!$OMP PARALLEL
143!$OMP WORKSHARE
144! calcul de bdx et hdx
145hdx(:,:)=dx1*(vieuxH(:,:)-eoshift(vieuxH(:,:),shift=-1,boundary=0d0,dim=1))
146hdy(:,:)=dx1*(vieuxH(:,:)-eoshift(vieuxH(:,:),shift=-1,boundary=0d0,dim=2))
147bdx(:,:)=dx1*(B(:,:)-eoshift(B(:,:),shift=-1,boundary=0d0,dim=1))
148bdy(:,:)=dx1*(B(:,:)-eoshift(B(:,:),shift=-1,boundary=0d0,dim=2))
149
150
151! initialisations (qui feront aussi les conditions aux limites)
152Arelax(:,:)=0.
153Brelax(:,:)=0.
154Drelax(:,:)=0.
155Erelax(:,:)=0.
156Crelax(:,:)=1.
157Frelax(:,:)=0.
158DeltaH(:,:)=0.
159!$OMP END WORKSHARE
160
161! schema spatial
162
163if (upwind.eq.1.) then                 !schema amont
164!$OMP WORKSHARE
165   where (Advx(:,:).ge.0.)
166      c_west(:,:)=1.
167      c_east(:,:)=0.
168   elsewhere
169      c_west(:,:)=0.
170      c_east(:,:)=1.
171   end where
172
173   where (Advy(:,:).ge.0.)
174      c_south(:,:)=1.
175      c_north(:,:)=0.
176   elsewhere
177      c_south(:,:)=0.
178      c_north(:,:)=1.
179   end where
180!$OMP END WORKSHARE
181else if (upwind.lt.1.) then             ! schema centre
182!$OMP WORKSHARE
183      c_west(:,:)=0.5
184      c_east(:,:)=0.5
185      c_south(:,:)=0.5
186      c_north(:,:)=0.5
187!$OMP END WORKSHARE     
188end if
189
190
191! attribution des elements des diagonales
192!$OMP DO PRIVATE(frdx,frdy,fraxw,fraxe,frays,frayn)
193do j=2,ny-1
194  do i=2,nx-1
195
196!  sous diagonale en x
197     arelax(i,j)=-omega*Dtdx2*Dfx(i,j)   &               ! partie diffusive en x
198          -mu_adv*dtdx*c_west(i,j)*Advx(i,j)             ! partie advective en x
199
200!  sur diagonale en x
201     brelax(i,j)=-omega*Dtdx2*Dfx(i+1,j)  &              ! partie diffusive
202          +mu_adv*dtdx*c_east(i+1,j)*Advx(i+1,j)         ! partie advective
203
204!  sous diagonale en y
205     drelax(i,j)=-omega*Dtdx2*Dfy(i,j)   &               ! partie diffusive en y
206          -mu_adv*dtdx*c_south(i,j)*Advy(i,j)            ! partie advective en y
207
208!  sur diagonale en y
209     erelax(i,j)=-omega*Dtdx2*Dfy(i,j+1)  &              ! partie diffusive
210          +mu_adv*dtdx*c_north(i,j+1)*Advy(i,j+1)        ! partie advective
211
212
213
214! diagonale
215     crelax(i,j)=omega*Dtdx2*((Dfx(i+1,j)+Dfx(i,j))  &                         ! partie diffusive en x
216                             +(Dfy(i,j+1)+Dfy(i,j)))                           ! partie diffusive en y
217     crelax(i,j)=crelax(i,j)+mu_adv*dtdx* &
218                  ( (c_west(i+1,j)*Advx(i+1,j)-c_east(i,j)*Advx(i,j)) &        !partie advective en x
219                  +(c_south(i,j+1)*Advy(i,j+1)-c_north(i,j)*Advy(i,j)))        !partie advective en y
220     crelax(i,j)=1.+crelax(i,j)                                                ! partie temporelle
221
222
223! terme du vecteur
224
225     frdx= -Dfx(i,j) * (Bdx(i,j)  +(1.-omega)*Hdx(i,j))         &  ! partie diffusive en x
226           +Dfx(i+1,j)*(Bdx(i+1,j)+(1.-omega)*Hdx(i+1,j))
227
228     frdy= -Dfy(i,j) * (Bdy(i,j)  +(1.-omega)*Hdy(i,j))         &  ! partie diffusive en y
229           +Dfy(i,j+1)*(Bdy(i,j+1)+(1.-omega)*Hdy(i,j+1))
230
231     fraxw= -c_west(i,j)*  Advx(i,j) * vieuxH(i-1,j)           &  ! partie advective en x
232            +c_west(i+1,j)*Advx(i+1,j)*vieuxH(i,j)                ! venant de l'west
233
234     fraxe= -c_east(i,j) * Advx(i,j) * vieuxH(i,j)             &  ! partie advective en x
235            +c_east(i+1,j)*Advx(i+1,j)*vieuxH(i+1,j)              ! venant de l'est
236
237     frays= -c_south(i,j) * Advy(i,j) * vieuxH(i,j-1)          &  ! partie advective en y
238            +c_south(i,j+1)*Advy(i,j+1)*vieuxH(i,j)               ! venant du sud
239
240     frayn= -c_north(i,j) * Advy(i,j) * vieuxH(i,j)            &  ! partie advective en y
241            +c_north(i,j+1)*Advy(i,j+1)*vieuxH(i,j+1)             ! venant du nord
242
243
244
245     frelax(i,j) = vieuxH(i,j) + Dt*bil(i,j)    &
246                 + dtdx*(frdx+frdy)        &  ! termes lies a la diffusion (socle, omega)
247                 + (1.-mu_adv)*dtdx*((fraxw+fraxe)+(frays+frayn)) 
248                                              ! le dernier terme serait pour
249                                              ! de l'over implicite en diffusion
250                                              ! (en pratique pas utilise)
251     end do
252  end do
253!$OMP END DO
254
255!$OMP WORKSHARE
256where (i_hpresc(:,:) .eq.1)          ! thickness prescribed
257   frelax(:,:) = H_presc(:,:)
258   arelax(:,:) = 0.
259   brelax(:,:) = 0.
260   crelax(:,:) = 1.
261   drelax(:,:) = 0.
262   erelax(:,:) = 0.
263end where
264!$OMP END WORKSHARE
265!$OMP END PARALLEL
266
267stopp = .false.
268ntour=0
269
270relax_loop: do while(.not.stopp)
271   ntour=ntour+1
272   !$OMP PARALLEL
273   !$OMP DO PRIVATE(reste)
274   do j=2,ny-1
275      do i=2,nx-1
276         reste = (((arelax(i,j)*newH(i-1,j) +drelax(i,j)*newH(i,j-1)) &
277               + (brelax(i,j)*newH(i+1,j) + erelax(i,j)*newH(i,j+1))) &
278               + crelax(i,j)*newH(i,j))- frelax(i,j)
279
280         deltaH(i,j) = reste/crelax(i,j)
281      end do
282   end do
283   !$OMP END DO
284
285   !$OMP WORKSHARE
286   newH(:,:)=newH(:,:)-deltaH(:,:)
287   !$OMP END WORKSHARE
288   !$OMP END PARALLEL
289
290
291! critere d'arret:
292! ----------------         
293   delh=0
294
295   !$OMP PARALLEL
296   !$OMP DO REDUCTION(+:delh)
297   do j=2,ny-1
298      do i=2,nx-1
299         delh=delh+deltaH(i,j)**2
300      end do
301   end do
302   !$OMP END DO
303   !$OMP END PARALLEL
304
305   if (delh.gt.0.) then
306      testh=sqrt(delh)/((nx-2)*(ny-2))
307   else
308      testh=0.
309   endif
310   stopp = (testh.lt.1.e-4).or.(ntour.gt.100)
311
312end do relax_loop
313
314if (itracebug.eq.1)  call tracebug(' fin routine resolution_diffusion')
315 
316return
317
318end subroutine resol_adv_diff_2D_vect
319
320end module reso_adv_diff_2D_vect
321
322
323
324
Note: See TracBrowser for help on using the repository browser.