source: branches/iLoveclim/SOURCES/resol_adv_diff_2D-sept2009.f90 @ 146

Last change on this file since 146 was 146, checked in by aquiquet, 7 years ago

Grisli-iLoveclim branch: merged to trunk at revision 145

File size: 10.6 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
103!$ USE OMP_LIB
104implicit none
105
106real,dimension(nx,ny), intent(in) :: Dfx          !< terme diffusif selon x
107real,dimension(nx,ny), intent(in) :: Dfy          !< terme diffusif selon y
108real,dimension(nx,ny), intent(in) :: Advx         !< terme advectif selon x
109real,dimension(nx,ny), intent(in) :: Advy         !< terme advectif selon y
110real,dimension(nx,ny), intent(in) :: vieuxH       !< ancienne valeur de H
111real,dimension(nx,ny), intent(out):: newH         !< nouvelle valeur de H
112
113real,dimension(nx,ny), intent(in) :: bil          !< bilan de masse pour la colonne
114
115real,dimension(nx,ny), intent(in) :: H_presc      !< H value if prescribed
116integer,dimension(nx,ny), intent(in) :: i_Hpresc  !< 1 if H is prescribed on this node, else 0
117
118
119! tableaux de travail. resolution M H = Frelax
120!!$real,dimension(nx,ny) :: crelax      !< diagnonale de M
121!!$real,dimension(nx,ny) :: arelax      !< sous diagonale selon x
122!!$real,dimension(nx,ny) :: brelax      !< sur diagonale selon x
123!!$real,dimension(nx,ny) :: drelax      !< sous diagonale selon y
124!!$real,dimension(nx,ny) :: erelax      !< sur diagonale selon y
125!!$real,dimension(nx,ny) :: frelax      !< vecteur
126!!$real,dimension(nx,ny) :: c_west      !< sur demi mailles Ux
127!!$real,dimension(nx,ny) :: c_east      !< sur demi mailles Ux
128!!$real,dimension(nx,ny) :: c_north     !< sur demi mailles Uy
129!!$real,dimension(nx,ny) :: c_south     !<sur demi mailles Ux
130
131!!$real,dimension(nx,ny) :: bdx         !< pente socle
132!!$real,dimension(nx,ny) :: bdy         !< pente socle
133!!$
134!!$real,dimension(nx,ny) :: hdx         !< pente epaisseur
135!!$real,dimension(nx,ny) :: hdy         !< pente epaisseur
136
137real :: frdx,frdy                    !< pour calcul frelax : termes diffusion
138real :: fraxw,fraxe,frays,frayn      !< termes advection
139
140real,dimension(nx,ny) :: deltah      ! dans calcul relax
141real :: delh                         ! dans calcul relax
142real :: testh                        ! dans calcul relax
143
144logical :: stopp
145integer :: ntour
146real :: reste
147
148
149if (itracebug.eq.1)  call tracebug(' Entree dans routine resolution_diffusion')
150
151!$OMP PARALLEL
152!$OMP WORKSHARE
153! calcul de bdx et hdx
154hdx(:,:)=dx1*(vieuxH(:,:)-eoshift(vieuxH(:,:),shift=-1,boundary=0.,dim=1))
155hdy(:,:)=dx1*(vieuxH(:,:)-eoshift(vieuxH(:,:),shift=-1,boundary=0.,dim=2))
156bdx(:,:)=dx1*(B(:,:)-eoshift(B(:,:),shift=-1,boundary=0.,dim=1))
157bdy(:,:)=dx1*(B(:,:)-eoshift(B(:,:),shift=-1,boundary=0.,dim=2))
158
159
160! initialisations (qui feront aussi les conditions aux limites)
161Arelax(:,:)=0.
162Brelax(:,:)=0.
163Drelax(:,:)=0.
164Erelax(:,:)=0.
165Crelax(:,:)=1.
166Frelax(:,:)=0.
167DeltaH(:,:)=0.
168!$OMP END WORKSHARE
169
170! schema spatial
171
172if (upwind.eq.1.) then                 !schema amont
173!$OMP WORKSHARE
174   where (Advx(:,:).ge.0.)
175      c_west(:,:)=1.
176      c_east(:,:)=0.
177   elsewhere
178      c_west(:,:)=0.
179      c_east(:,:)=1.
180   end where
181
182   where (Advy(:,:).ge.0.)
183      c_south(:,:)=1.
184      c_north(:,:)=0.
185   elsewhere
186      c_south(:,:)=0.
187      c_north(:,:)=1.
188   end where
189!$OMP END WORKSHARE
190else if (upwind.lt.1.) then             ! schema centre
191!$OMP WORKSHARE
192      c_west(:,:)=0.5
193      c_east(:,:)=0.5
194      c_south(:,:)=0.5
195      c_north(:,:)=0.5
196!$OMP END WORKSHARE     
197end if
198
199
200! attribution des elements des diagonales
201!$OMP DO PRIVATE(frdx,frdy,fraxw,fraxe,frays,frayn)
202do j=2,ny-1
203  do i=2,nx-1
204
205!  sous diagonale en x
206     arelax(i,j)=-omega*Dtdx2*Dfx(i,j)   &               ! partie diffusive en x
207          -mu_adv*dtdx*c_west(i,j)*Advx(i,j)             ! partie advective en x
208
209!  sur diagonale en x
210     brelax(i,j)=-omega*Dtdx2*Dfx(i+1,j)  &              ! partie diffusive
211          +mu_adv*dtdx*c_east(i+1,j)*Advx(i+1,j)         ! partie advective
212
213!  sous diagonale en y
214     drelax(i,j)=-omega*Dtdx2*Dfy(i,j)   &               ! partie diffusive en y
215          -mu_adv*dtdx*c_south(i,j)*Advy(i,j)            ! partie advective en y
216
217!  sur diagonale en y
218     erelax(i,j)=-omega*Dtdx2*Dfy(i,j+1)  &              ! partie diffusive
219          +mu_adv*dtdx*c_north(i,j+1)*Advy(i,j+1)        ! partie advective
220
221
222
223! diagonale
224     crelax(i,j)=omega*Dtdx2*((Dfx(i+1,j)+Dfx(i,j))  &                         ! partie diffusive en x
225                             +(Dfy(i,j+1)+Dfy(i,j)))                           ! partie diffusive en y
226     crelax(i,j)=crelax(i,j)+mu_adv*dtdx* &
227                  ( (c_west(i+1,j)*Advx(i+1,j)-c_east(i,j)*Advx(i,j)) &        !partie advective en x
228                  +(c_south(i,j+1)*Advy(i,j+1)-c_north(i,j)*Advy(i,j)))        !partie advective en y
229     crelax(i,j)=1.+crelax(i,j)                                                ! partie temporelle
230
231
232! terme du vecteur
233
234     frdx= -Dfx(i,j) * (Bdx(i,j)  +(1.-omega)*Hdx(i,j))         &  ! partie diffusive en x
235           +Dfx(i+1,j)*(Bdx(i+1,j)+(1.-omega)*Hdx(i+1,j))
236
237     frdy= -Dfy(i,j) * (Bdy(i,j)  +(1.-omega)*Hdy(i,j))         &  ! partie diffusive en y
238           +Dfy(i,j+1)*(Bdy(i,j+1)+(1.-omega)*Hdy(i,j+1))
239
240     fraxw= -c_west(i,j)*  Advx(i,j) * vieuxH(i-1,j)           &  ! partie advective en x
241            +c_west(i+1,j)*Advx(i+1,j)*vieuxH(i,j)                ! venant de l'west
242
243     fraxe= -c_east(i,j) * Advx(i,j) * vieuxH(i,j)             &  ! partie advective en x
244            +c_east(i+1,j)*Advx(i+1,j)*vieuxH(i+1,j)              ! venant de l'est
245
246     frays= -c_south(i,j) * Advy(i,j) * vieuxH(i,j-1)          &  ! partie advective en y
247            +c_south(i,j+1)*Advy(i,j+1)*vieuxH(i,j)               ! venant du sud
248
249     frayn= -c_north(i,j) * Advy(i,j) * vieuxH(i,j)            &  ! partie advective en y
250            +c_north(i,j+1)*Advy(i,j+1)*vieuxH(i,j+1)             ! venant du nord
251
252
253
254     frelax(i,j) = vieuxH(i,j) + Dt*bil(i,j)    &
255                 + dtdx*(frdx+frdy)        &  ! termes lies a la diffusion (socle, omega)
256                 + (1.-mu_adv)*dtdx*((fraxw+fraxe)+(frays+frayn)) 
257                                              ! le dernier terme serait pour
258                                              ! de l'over implicite en diffusion
259                                              ! (en pratique pas utilise)
260     end do
261  end do
262!$OMP END DO
263
264!$OMP WORKSHARE
265where (i_hpresc(:,:) .eq.1)          ! thickness prescribed
266   frelax(:,:) = H_presc(:,:)
267   arelax(:,:) = 0.
268   brelax(:,:) = 0.
269   crelax(:,:) = 1.
270   drelax(:,:) = 0.
271   erelax(:,:) = 0.
272end where
273!$OMP END WORKSHARE
274!$OMP END PARALLEL
275
276stopp = .false.
277ntour=0
278
279relax_loop: do while(.not.stopp)
280   ntour=ntour+1
281   !$OMP PARALLEL
282   !$OMP DO PRIVATE(reste)
283   do j=2,ny-1
284      do i=2,nx-1
285         reste = (((arelax(i,j)*newH(i-1,j) +drelax(i,j)*newH(i,j-1)) &
286               + (brelax(i,j)*newH(i+1,j) + erelax(i,j)*newH(i,j+1))) &
287               + crelax(i,j)*newH(i,j))- frelax(i,j)
288
289         deltaH(i,j) = reste/crelax(i,j)
290      end do
291   end do
292   !$OMP END DO
293
294   !$OMP WORKSHARE
295   newH(:,:)=newH(:,:)-deltaH(:,:)
296   !$OMP END WORKSHARE
297   !$OMP END PARALLEL
298
299
300! critere d'arret:
301! ----------------         
302   delh=0
303
304   !$OMP PARALLEL
305   !$OMP DO REDUCTION(+:delh)
306   do j=2,ny-1
307      do i=2,nx-1
308         delh=delh+deltaH(i,j)**2
309      end do
310   end do
311   !$OMP END DO
312   !$OMP END PARALLEL
313
314   if (delh.gt.0.) then
315      testh=sqrt(delh)/((nx-2)*(ny-2))
316   else
317      testh=0.
318   endif
319   stopp = (testh.lt.1.e-4).or.(ntour.gt.100)
320
321end do relax_loop
322
323if (itracebug.eq.1)  call tracebug(' fin routine resolution_diffusion')
324 
325return
326
327end subroutine resol_adv_diff_2D_vect
328
329end module reso_adv_diff_2D_vect
330
331
332
333
Note: See TracBrowser for help on using the repository browser.