source: trunk/LMDZ.COMMON/libf/evolution/math_mod.F90 @ 3599

Last change on this file since 3599 was 3553, checked in by jbclement, 6 weeks ago

PEM:
Addition of the shifting of the soil temperature profile to follow the surface evolution due to ice condensation/sublimation + small cleanings/improvements.
JBC

File size: 8.8 KB
RevLine 
[2961]1module math_mod
2!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3!!!
[3526]4!!! Purpose: The module contains all the mathematical SUBROUTINE used in the PEM
[3532]5!!!
[2961]6!!! Author: Adapted from Schorgofer MSIM (N.S, Icarus 2010), impletented here by LL
7!!!
8!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[3526]9
10implicit none
11
12!=======================================================================
[2961]13  contains
[3526]14!=======================================================================
[2961]15
[3526]16SUBROUTINE deriv1(z,nz,y,y0,ybot,dzY)
[2961]17!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
18!!!
[3532]19!!! Purpose: Compute the first derivative of a function y(z) on an irregular grid
20!!!
[2961]21!!! Author: From N.S (N.S, Icarus 2010), impletented here by LL
22!!!
23!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
24! first derivative of a function y(z) on irregular grid
[3526]25! upper boundary conditions: y(0) = y0
[2961]26! lower boundary condition.: yp = ybottom
[3526]27
28implicit none
29
[3553]30! Inputs
31!-------
32integer,             intent(in) :: nz       ! number of layer
33real, dimension(nz), intent(in) :: z        ! depth layer
34real, dimension(nz), intent(in) :: y        ! function which needs to be derived
35real,                intent(in) :: y0, ybot ! boundary conditions
36! Outputs
37!--------
38real, dimension(nz), intent(out) :: dzY ! derivative of y w.r.t depth
39! Local variables
40!----------------
[3526]41integer :: j
42real    :: hm, hp, c1, c2, c3
[2961]43
[3526]44hp = z(2) - z(1)
45c1 = z(1)/(hp*z(2))
46c2 = 1/z(1) - 1/(z(2) - z(1))
47c3 = -hp/(z(1)*z(2))
48dzY(1) = c1*y(2) + c2*y(1) + c3*y0
49do j = 2,nz - 1
50    hp = z(j + 1) - z(j)
51    hm = z(j) - z(j - 1)
52    c1 = +hm/(hp*(z(j + 1) - z(j - 1)))
53    c2 = 1/hm - 1/hp
54    c3 = -hp/(hm*(z(j + 1) - z(j - 1)))
55    dzY(j) = c1*y(j + 1) + c2*y(j) + c3*y(j - 1)
56enddo
57dzY(nz) = (ybot - y(nz - 1))/(2.*(z(nz) - z(nz - 1)))
[2961]58
[3526]59END SUBROUTINE deriv1
[2961]60
[3526]61!=======================================================================
[2961]62
[3526]63SUBROUTINE deriv2_simple(z,nz,y,y0,yNp1,yp2)
[2961]64!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
65!!!
[3532]66!!! Purpose: Compute the second derivative of a function y(z) on an irregular grid
67!!!
[2961]68!!! Author: N.S (raw copy/paste from MSIM)
69!!!
70!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[3526]71! second derivative y_zz on irregular grid
72! boundary conditions: y(0) = y0, y(nz + 1) = yNp1
[2961]73
[3526]74implicit none
[2961]75
[3553]76! Inputs
77!-------
78integer,             intent(in) :: nz
79real, dimension(nz), intent(in) :: z, y
80real,                intent(in) :: y0, yNp1
81! Outputs
82!--------
83real, dimension(nz), intent(out) :: yp2
84! Local variables
85!----------------
[3526]86integer :: j
87real    :: hm, hp, c1, c2, c3
[2961]88
[3526]89c1 = +2./((z(2) - z(1))*z(2))
90c2 = -2./((z(2) - z(1))*z(1))
91c3 = +2./(z(1)*z(2))
92yp2(1) = c1*y(2) + c2*y(1) + c3*y0
93do j = 2,nz - 1
94    hp = z(j + 1) - z(j)
95    hm = z(j) - z(j - 1)
96    c1 = +2./(hp*(z(j + 1) - z(j - 1)))
97    c2 = -2./(hp*hm)
98    c3 = +2./(hm*(z(j + 1) - z(j-1)))
99    yp2(j) = c1*y(j + 1) + c2*y(j) + c3*y(j - 1)
100enddo
101yp2(nz) = (yNp1 - 2*y(nz) + y(nz - 1))/(z(nz) - z(nz - 1))**2
[2961]102
[3526]103END SUBROUTINE deriv2_simple
[2961]104
[3526]105!=======================================================================
[2961]106
[3526]107SUBROUTINE  deriv1_onesided(j,z,nz,y,dy_zj)
[2961]108!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
109!!!
[3526]110!!! Purpose: First derivative of function y(z) at z(j)  one-sided derivative on irregular grid
[3532]111!!!
[2961]112!!! Author: N.S (raw copy/paste from MSIM)
113!!!
114!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
115
[3553]116implicit none
[2961]117
[3553]118! Inputs
119!-------
120integer,             intent(in) :: nz, j
121real, dimension(nz), intent(in) :: z, y
122! Outputs
123!--------
[3526]124real, intent(out) :: dy_zj
[3553]125! Local viariables
126!-----------------
[3526]127real :: h1, h2, c1, c2, c3
[2961]128
[3526]129if (j < 1 .or. j > nz - 2) then
130    dy_zj = -1.
131else
132    h1 = z(j + 1) - z(j)
133    h2 = z(j + 2)- z(j + 1)
134    c1 = -(2*h1 + h2)/(h1*(h1 + h2))
135    c2 = (h1 + h2)/(h1*h2)
136    c3 = -h1/(h2*(h1 + h2))
137    dy_zj = c1*y(j) + c2*y(j + 1) + c3*y(j + 2)
138endif
139
140END SUBROUTINE deriv1_onesided
141
142!=======================================================================
143
144PURE SUBROUTINE colint(y,z,nz,i1,i2,integral)
[2961]145!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
146!!!
147!!! Purpose:  Column integrates y on irregular grid
[3532]148!!!
[2961]149!!! Author: N.S (raw copy/paste from MSIM)
150!!!
151!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[3532]152
[3526]153implicit none
[3532]154
[3553]155! Inputs
156!-------
157integer,             intent(in) :: nz, i1, i2
158real, dimension(nz), intent(in) :: y, z
159! Outputs
160!--------
[3526]161real, intent(out) :: integral
[3553]162! Local viariables
163!-----------------
164integer             :: i
165real, dimension(nz) :: dz
[2961]166
[3526]167dz(1) = (z(2) - 0.)/2
168do i = 2,nz - 1
169    dz(i) = (z(i + 1) - z(i - 1))/2.
170enddo
171dz(nz) = z(nz) - z(nz - 1)
172integral = sum(y(i1:i2)*dz(i1:i2))
[2961]173
[3526]174END SUBROUTINE colint
[2961]175
[3526]176!=======================================================================
[2961]177
178SUBROUTINE findroot(y1,y2,z1,z2,zr)
[3526]179!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
180!!!
181!!! Purpose: Compute the root zr, between two values y1 and y2 at depth z1,z2
[3532]182!!!
[3526]183!!! Author: LL
184!!!
185!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[2961]186
[3526]187implicit none
[2961]188
[3553]189! Inputs
190!-------
[3526]191real, intent(in) :: y1, y2 ! difference between surface water density and at depth  [kg/m^3]
192real, intent(in) :: z1, z2 ! depth [m]
[3553]193! Outputs
194!--------
195real, intent(out) :: zr ! depth at which we have zero
[3526]196
197zr = (y1*z2 - y2*z1)/(y1 - y2)
198
199END SUBROUTINE findroot
200
[3553]201!=======================================================================
202
203SUBROUTINE solve_tridiag(a,b,c,d,n,x,error)
204!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
205!!!
206!!! Purpose: Solve a tridiagonal system Ax = d using the Thomas' algorithm
207!!!          a: sub-diagonal
208!!!          b: main diagonal
209!!!          c: super-diagonal
210!!!          d: right-hand side
211!!!          x: solution
212!!!
213!!! Author: JBC
214!!!
215!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
216
217implicit none
218
219! Inputs
220!-------
221integer,                intent(in) :: n
222real, dimension(n),     intent(in) :: b, d
223real, dimension(n - 1), intent(in) :: a, c
224! Outputs
225!--------
226real, dimension(n), intent(out) :: x
227integer,            intent(out) :: error
228! Local viariables
229!-----------------
230integer            :: i
231real               :: m
232real, dimension(n) :: cp, dp
233
234! Check stability: diagonally dominant condition
235error = 0
236if (abs(b(1)) < abs(c(1))) then
237    error = 1
238    return
239endif
240do i = 2,n - 1
241    if (abs(b(i)) < abs(a(i - 1)) + abs(c(i))) then
242        error = 1
243        return
244    endif
245enddo
246if (abs(b(n)) < abs(a(n - 1))) then
247    error = 1
248    return
249endif
250
251! Initialization
252cp(1) = c(1)/b(1)
253dp(1) = d(1)/b(1)
254
255! Forward sweep
256do i = 2,n - 1
257    m = b(i) - a(i - 1)*cp(i - 1)
258    cp(i) = c(i)/m
259    dp(i) = (d(i) - a(i - 1)*dp(i - 1))/m
260enddo
261m = b(n) - a(n - 1)*cp(n - 1)
262dp(n) = (d(n) - a(n - 1)*dp(n - 1))/m
263
264! Backward substitution
265x(n) = dp(n)
266do i = n - 1,1,-1
267    x(i) = dp(i) - cp(i)*x(i + 1)
268enddo
269
270END SUBROUTINE solve_tridiag
271
272!=======================================================================
273
274SUBROUTINE solve_steady_heat(n,z,mz,kappa,mkappa,T_left,q_right,T)
275!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
276!!!
277!!! Purpose: Solve 1D steady-state heat equation with space-dependent diffusivity
278!!!          Left boudary condition : prescribed temperature T_left
279!!!          Right boudary condition: prescribed thermal flux q_right
280!!!
281!!!          z     : grid points
282!!!          mz    : mid-grid points
283!!!          kappa : thermal diffusivity at grid points
284!!!          mkappa: thermal diffusivity at mid-grid points
285!!!
286!!! Author: JBC
287!!!
288!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
289
290use abort_pem_mod, only: abort_pem
291
292implicit none
293
294! Inputs
295!-------
296integer,                intent(in) :: n
297real, dimension(n),     intent(in) :: z, kappa
298real, dimension(n - 1), intent(in) :: mz, mkappa
299real,                   intent(in) :: T_left, q_right
300! Outputs
301!--------
302real, dimension(n), intent(out) :: T
303! Local viariables
304!-----------------
305integer                :: i, error
306real, dimension(n)     :: b, d
307real, dimension(n - 1) :: a, c
308
309! Initialization
310a = 0.; b = 0.; c = 0.; d = 0.
311
312! Left boundary condition (Dirichlet: prescribed temperature)
313b(1) = 1.
314d(1) = T_left
315
316! Internal points
317do i = 2,n - 1
318    a(i - 1) = -mkappa(i - 1)/((mz(i) - mz(i - 1))*(z(i) - z(i - 1)))
319    c(i) = -mkappa(i)/((mz(i) - mz(i - 1))*(z(i + 1) - z(i)))
320    b(i) = -(a(i - 1) + c(i))
321enddo
322
323! Right boundary condition (Neumann: prescribed temperature)
324a(n - 1) = kappa(n - 1)/(z(n) - z(n - 1))
325b(n) = -kappa(n)/(z(n) - z(n - 1))
326d(n) = q_right
327
328! Solve the tridiagonal linear system with the Thomas' algorithm
329call solve_tridiag(a,b,c,d,n,T,error)
330if (error /= 0) call abort_pem("solve_steady_heat","Unstable solving!",1)
331
332END SUBROUTINE solve_steady_heat
333
[2961]334end module math_mod
Note: See TracBrowser for help on using the repository browser.