source: trunk/SOURCES/relaxation_water_diffusion.f90 @ 111

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

OpenMP parallelization in conserv-mass-adv-diff_sept2009_mod.f90, diffusiv-polyn-0.6.f90, dragging_neff_slope_mod.f90, eaubasale-0.5_mod.f90 and relaxation_water_diffusion.f90

File size: 6.5 KB
Line 
1!> \file relaxation_water_diffusion.f90
2!! Module pour la resolution de l'equation de relaxation et de l equation de diffusion.
3!<
4
5!> \namespace relaxation_waterdif_mod
6!! Module pour la resolution de l'equation de relaxation et de l equation de diffusion .
7!! \author ...
8!! \date ...
9!! @note dhwat/dt = bmelt-infiltr-d/dx(Kond*dhwat/dx)+d/dx(Kond*pgx)
10!<
11module relaxation_waterdif_mod
12
13
14
15
16CONTAINS
17
18
19subroutine relaxation_waterdif(NXX,NYY,DT,DX,vieuxHWATER,limit_hw,klimit,BMELT,INFILTR,PGMX,PGMY,KOND,KONDMAX,HWATER)
20
21  !$ USE OMP_LIB
22
23  implicit none
24
25
26    ! declaration des variables en entree
27    !------------------------------------------------
28    INTEGER, intent(in) :: NXX, NYY       !< defini la taille des tableaux
29    REAL,    intent(in) ::  DT            !< pas de temps court
30    REAL,    intent(in) ::  DX            !< pas en x
31    REAL,    intent(in) :: INFILTR        !< basal infiltration (lose of water)
32    REAL,    intent(in) :: KONDMAX        !< maximum hydaulic conductivity (outside ice sheet)
33
34    REAL,dimension(NXX,NYY), intent(in) :: limit_hw    !< conditions aux limites
35    integer,dimension(NXX,NYY), intent(in) :: klimit    !< ou appliquer les conditions
36    REAL,dimension(NXX,NYY), intent(in) :: vieuxHWATER    !< H au pas de temps precedent 'o'
37    REAL,dimension(NXX,NYY), intent(in) :: BMELT     !< basal melting  'o'
38    REAL,dimension(NXX,NYY), intent(in) :: PGMX     !< hydaulic potential gratient '>'
39    REAL,dimension(NXX,NYY), intent(in) :: PGMY     !< hydaulic potential gratient '^'
40    REAL,dimension(NXX,NYY), intent(in) :: KOND     !< hydaulic conductivity 'o'
41
42    ! declaration des variables en sortie
43    !-------------------------------------
44    REAL,dimension(NXX,NYY), intent(out):: HWATER      !< basal water thickness  'o' (pressure equivalent)
45
46
47    ! declaration des variables locales
48    !----------------------------------
49    INTEGER :: I,J
50    REAL  :: TESTH 
51    REAL,dimension(NXX,NYY) :: ARELAX,BRELAX,CRELAX,DRELAX,ERELAX,FRELAX
52    REAL,dimension(NXX,NYY) :: DELTAH
53    REAL :: RESTE,DELH,VH
54    INTEGER :: ntour
55    REAL  :: DTSRGDX,dtwdx2
56    LOGICAL :: STOPP
57    REAL,dimension(NXX,NYY) :: KMX, KMY
58    !    REAL  :: RHO,RHOW,RHOG  !,SECYEAR!! ice density, water density, density*acceleration
59!    write(6,*) 'entree relaxation'
60
61! write(166,*)' entree  relaxation waterdif'
62!$OMP PARALLEL
63!$OMP WORKSHARE
64    HWATER(:,:)= vieuxHWATER(:,:)
65!$OMP END WORKSHARE
66    !   calcul de kmx et kmx a partir de KOND
67    !   conductivite hyrdraulique sur les noeuds mineurs
68    !   moyenne harmonique
69    !   ----------------------------------------
70!$OMP DO
71    do j=2,nyy
72       do i=2,nxx
73
74          if  ((kond(i,j).lt.1.e-20).or.(kond(i-1,j).lt.1.e-20)) then
75             kmx(i,j)=0. ! to avoid division by o
76          else
77             kmx(i,j)=2*(kond(i,j)*kond(i-1,j))/(kond(i,j)+kond(i-1,j))
78          endif
79
80       end do
81    end do
82!$OMP END DO
83
84!$OMP DO
85    do j=2,nyy
86       do i=2,nxx
87          if  ((kond(i,j).lt.1.e-20).or.(kond(i,j-1).lt.1.e-20)) then
88             kmy(i,j)=0. ! to avoid division by o
89          else
90             kmy(i,j)=2*(kond(i,j)*kond(i,j-1))/(kond(i,j)+kond(i,j-1))
91          endif
92
93       enddo
94    enddo
95!$OMP END DO
96
97    !   attribution des coefficients  arelax ....
98    !   ----------------------------------------
99    !    SECYEAR=365.*24.*3600.
100    !   rho=910.
101    !    rhow=1000.
102    !    rhog=rhow*9.81
103    !    dtsrgdx=dt/(rhog*DX) a mon avis c'est rhow qu'il fallait utiliser. Maintenant cette
104    !    division est faite dans eaubasale
105
106    dtsrgdx=dt/DX
107    dtwdx2=dt/dx/dx
108
109!$OMP WORKSHARE
110  arelax(:,:)=0.
111  brelax(:,:)=0.
112  crelax(:,:)=1.
113  drelax(:,:)=0.
114  erelax(:,:)=0.
115  frelax(:,:)=limit_hw(:,:)
116!$OMP END WORKSHARE
117
118!$OMP DO
119    do J=2,NYY-1
120       do I=2,NXX-1
121
122          if (klimit(i,j).eq.0) then
123! calcul du vecteur
124          FRELAX(I,J)= VIEUXHWATER(I,J)+(BMELT(I,J)-INFILTR)*DT
125          frelax(i,j)=frelax(i,j)+(kmx(i,j)*pgmx(i,j)-kmx(i+1,j)*pgmx(i+1,j))*dtsrgdx
126          frelax(i,j)=frelax(i,j)+(kmy(i,j)*pgmy(i,j)-kmy(i,j+1)*pgmy(i,j+1))*dtsrgdx
127! calcul des diagonales     
128          arelax(i,j)=-kmx(i,j)*dtwdx2        ! arelax : diagonale i-1,j
129
130          brelax(i,j)=-kmx(i+1,j)*dtwdx2      ! brelax : diagonale i+1,j
131
132          drelax(i,j)=-kmy(i,j)*dtwdx2        ! drelax : diagonale i,j-1
133         
134          erelax(i,j)=-kmy(i,j+1)*dtwdx2      ! drelax : diagonale i,j+1
135
136          crelax(i,j)=1.+((kmx(i,j)+kmx(i+1,j))+(kmy(i,j)+kmy(i,j+1)))*dtwdx2 
137                                              !crelax : diagonale i,j
138          else if (klimit(i,j).eq.1) then
139             hwater(i,j)=limit_hw(i,j)
140!             write(6,*) i,j,hwater(i,j),crelax(i,j),frelax(i,j),arelax(i,j)
141          endif
142       end do
143    end do
144!$OMP END DO
145!$OMP END PARALLEL
146
147    ! Boucle de relaxation :
148    ! ----------------------
149    STOPP = .false.
150
151    ntour=0
152    stopp=.false.
153    Do  WHILE(.NOT.STOPP)
154       ntour=ntour+1
155!       write(6,*) 'boucle de relaxation numero',ntour
156       !$OMP PARALLEL
157       !$OMP DO PRIVATE(reste)
158       do j=2,NYY-1
159          do i=2,NXX-1
160
161             RESTE = (((ARELAX(I,J)*HWATER(I-1,J) + BRELAX(I,J)*HWATER(I+1,J)) &
162                  + (DRELAX(I,J)*HWATER(I,J-1) + ERELAX(I,J)*HWATER(I,J+1))) &
163                  + CRELAX(I,J)*HWATER(I,J))- FRELAX(I,J)
164
165             DELTAH(I,J) = RESTE/CRELAX(I,J)             
166          end do
167       end do
168       !$OMP END DO
169
170       !$OMP WORKSHARE
171       deltah(:,:)=min(deltah(:,:),10.)
172       deltah(:,:)=max(deltah(:,:),-10.)
173       !$OMP END WORKSHARE
174       ! il faut faire le calcul suivant dans une autre boucle car RESTE est fonction
175       ! de hwater sur les points voisins.
176       !$OMP DO
177       do j=2,NYY-1
178          do i=2,NXX-1
179             HWATER(I,J) = HWATER(I,J) - DELTAH(I,J)
180          end do
181       end do
182       !$OMP END DO
183
184       ! critere d'arret:
185       Delh=0
186       Vh=0
187       
188       !$OMP DO REDUCTION(+:Delh)
189       DO j=2,NYY-1
190          DO i=2,NXX-1
191             
192             !   write(166,*) I,J,delh,deltah(i,j)
193             Delh=Delh+deltah(i,j)**2
194             !            Vh=Vh+h(i,j)**2.
195          END DO
196       END DO
197       !$OMP END DO
198       !$OMP END PARALLEL
199!       write(6,*) delh,maxval(deltah),minval(deltah)
200       !      testh=SQRT(Delh/Vh)
201       if (delh.gt.0.) then
202          testh=SQRT(Delh)/((NXX-2)*(NYY-2))
203       else
204          testh=0.
205       endif
206       STOPP = (testh.lt.1.E-3).or.(ntour.gt.1000)
207!       write(6,*) ntour,testh
208
209
210362   continue 
211end do
212  end subroutine relaxation_waterdif
213end module relaxation_waterdif_mod
Note: See TracBrowser for help on using the repository browser.