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