source: trunk/SOURCES/resol_adv_diff_2D-sept2009.f90 @ 65

Last change on this file since 65 was 65, checked in by dumas, 8 years ago

Deleting unused variables and move old sources

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