| 1 | module math_toolkit |
|---|
| 2 | !----------------------------------------------------------------------- |
|---|
| 3 | ! NAME |
|---|
| 4 | ! math_toolkit |
|---|
| 5 | ! |
|---|
| 6 | ! DESCRIPTION |
|---|
| 7 | ! The module contains all the mathematical subroutines used in the PEM. |
|---|
| 8 | ! |
|---|
| 9 | ! AUTHORS & DATE |
|---|
| 10 | ! L. Lange |
|---|
| 11 | ! JB Clement, 2023-2025 |
|---|
| 12 | ! |
|---|
| 13 | ! NOTES |
|---|
| 14 | ! Adapted from Schorghofer MSIM (N.S, Icarus 2010) |
|---|
| 15 | !----------------------------------------------------------------------- |
|---|
| 16 | |
|---|
| 17 | ! DECLARATION |
|---|
| 18 | ! ----------- |
|---|
| 19 | implicit none |
|---|
| 20 | |
|---|
| 21 | contains |
|---|
| 22 | !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|---|
| 23 | |
|---|
| 24 | !======================================================================= |
|---|
| 25 | SUBROUTINE deriv1(z,nz,y,y0,ybot,dzY) |
|---|
| 26 | !----------------------------------------------------------------------- |
|---|
| 27 | ! NAME |
|---|
| 28 | ! deriv1 |
|---|
| 29 | ! |
|---|
| 30 | ! DESCRIPTION |
|---|
| 31 | ! Compute the first derivative of a function y(z) on an irregular grid. |
|---|
| 32 | ! |
|---|
| 33 | ! AUTHORS & DATE |
|---|
| 34 | ! N.S (Icarus 2010) |
|---|
| 35 | ! L. Lange |
|---|
| 36 | ! JB Clement, 2023-2025 |
|---|
| 37 | ! |
|---|
| 38 | ! NOTES |
|---|
| 39 | ! Upper boundary conditions: y(0) = y0 |
|---|
| 40 | ! Lower boundary condition: yp = ybottom |
|---|
| 41 | !----------------------------------------------------------------------- |
|---|
| 42 | |
|---|
| 43 | ! DECLARATION |
|---|
| 44 | ! ----------- |
|---|
| 45 | implicit none |
|---|
| 46 | |
|---|
| 47 | ! ARGUMENTS |
|---|
| 48 | ! --------- |
|---|
| 49 | integer, intent(in) :: nz ! number of layer |
|---|
| 50 | real, dimension(nz), intent(in) :: z ! depth layer |
|---|
| 51 | real, dimension(nz), intent(in) :: y ! function which needs to be derived |
|---|
| 52 | real, intent(in) :: y0, ybot ! boundary conditions |
|---|
| 53 | real, dimension(nz), intent(out) :: dzY ! derivative of y w.r.t depth |
|---|
| 54 | |
|---|
| 55 | ! LOCAL VARIABLES |
|---|
| 56 | ! --------------- |
|---|
| 57 | integer :: j |
|---|
| 58 | real :: hm, hp, c1, c2, c3 |
|---|
| 59 | |
|---|
| 60 | ! CODE |
|---|
| 61 | ! ---- |
|---|
| 62 | hp = z(2) - z(1) |
|---|
| 63 | c1 = z(1)/(hp*z(2)) |
|---|
| 64 | c2 = 1/z(1) - 1/(z(2) - z(1)) |
|---|
| 65 | c3 = -hp/(z(1)*z(2)) |
|---|
| 66 | dzY(1) = c1*y(2) + c2*y(1) + c3*y0 |
|---|
| 67 | do j = 2,nz - 1 |
|---|
| 68 | hp = z(j + 1) - z(j) |
|---|
| 69 | hm = z(j) - z(j - 1) |
|---|
| 70 | c1 = +hm/(hp*(z(j + 1) - z(j - 1))) |
|---|
| 71 | c2 = 1/hm - 1/hp |
|---|
| 72 | c3 = -hp/(hm*(z(j + 1) - z(j - 1))) |
|---|
| 73 | dzY(j) = c1*y(j + 1) + c2*y(j) + c3*y(j - 1) |
|---|
| 74 | enddo |
|---|
| 75 | dzY(nz) = (ybot - y(nz - 1))/(2.*(z(nz) - z(nz - 1))) |
|---|
| 76 | |
|---|
| 77 | END SUBROUTINE deriv1 |
|---|
| 78 | !======================================================================= |
|---|
| 79 | |
|---|
| 80 | !======================================================================= |
|---|
| 81 | SUBROUTINE deriv2_simple(z,nz,y,y0,yNp1,yp2) |
|---|
| 82 | !----------------------------------------------------------------------- |
|---|
| 83 | ! NAME |
|---|
| 84 | ! deriv2_simple |
|---|
| 85 | ! |
|---|
| 86 | ! DESCRIPTION |
|---|
| 87 | ! Compute the second derivative of a function y(z) on an irregular grid. |
|---|
| 88 | ! |
|---|
| 89 | ! AUTHORS & DATE |
|---|
| 90 | ! N.S (Icarus 2010) |
|---|
| 91 | ! JB Clement, 2023-2025 |
|---|
| 92 | ! |
|---|
| 93 | ! NOTES |
|---|
| 94 | ! Boundary conditions: y(0) = y0, y(nz + 1) = yNp1 |
|---|
| 95 | !----------------------------------------------------------------------- |
|---|
| 96 | |
|---|
| 97 | ! DECLARATION |
|---|
| 98 | ! ----------- |
|---|
| 99 | implicit none |
|---|
| 100 | |
|---|
| 101 | ! ARGUMENTS |
|---|
| 102 | ! --------- |
|---|
| 103 | integer, intent(in) :: nz |
|---|
| 104 | real, dimension(nz), intent(in) :: z, y |
|---|
| 105 | real, intent(in) :: y0, yNp1 |
|---|
| 106 | real, dimension(nz), intent(out) :: yp2 |
|---|
| 107 | |
|---|
| 108 | ! LOCAL VARIABLES |
|---|
| 109 | ! --------------- |
|---|
| 110 | integer :: j |
|---|
| 111 | real :: hm, hp, c1, c2, c3 |
|---|
| 112 | |
|---|
| 113 | ! CODE |
|---|
| 114 | ! ---- |
|---|
| 115 | c1 = +2./((z(2) - z(1))*z(2)) |
|---|
| 116 | c2 = -2./((z(2) - z(1))*z(1)) |
|---|
| 117 | c3 = +2./(z(1)*z(2)) |
|---|
| 118 | yp2(1) = c1*y(2) + c2*y(1) + c3*y0 |
|---|
| 119 | do j = 2,nz - 1 |
|---|
| 120 | hp = z(j + 1) - z(j) |
|---|
| 121 | hm = z(j) - z(j - 1) |
|---|
| 122 | c1 = +2./(hp*(z(j + 1) - z(j - 1))) |
|---|
| 123 | c2 = -2./(hp*hm) |
|---|
| 124 | c3 = +2./(hm*(z(j + 1) - z(j-1))) |
|---|
| 125 | yp2(j) = c1*y(j + 1) + c2*y(j) + c3*y(j - 1) |
|---|
| 126 | enddo |
|---|
| 127 | yp2(nz) = (yNp1 - 2*y(nz) + y(nz - 1))/(z(nz) - z(nz - 1))**2 |
|---|
| 128 | |
|---|
| 129 | END SUBROUTINE deriv2_simple |
|---|
| 130 | !======================================================================= |
|---|
| 131 | |
|---|
| 132 | !======================================================================= |
|---|
| 133 | SUBROUTINE deriv1_onesided(j,z,nz,y,dy_zj) |
|---|
| 134 | !----------------------------------------------------------------------- |
|---|
| 135 | ! NAME |
|---|
| 136 | ! deriv1_onesided |
|---|
| 137 | ! |
|---|
| 138 | ! DESCRIPTION |
|---|
| 139 | ! First derivative of function y(z) at z(j) one-sided derivative |
|---|
| 140 | ! on irregular grid. |
|---|
| 141 | ! |
|---|
| 142 | ! AUTHORS & DATE |
|---|
| 143 | ! N.S (Icarus 2010) |
|---|
| 144 | ! JB Clement, 2023-2025 |
|---|
| 145 | ! |
|---|
| 146 | ! NOTES |
|---|
| 147 | ! |
|---|
| 148 | !----------------------------------------------------------------------- |
|---|
| 149 | |
|---|
| 150 | ! DECLARATION |
|---|
| 151 | ! ----------- |
|---|
| 152 | implicit none |
|---|
| 153 | |
|---|
| 154 | ! ARGUMENTS |
|---|
| 155 | ! --------- |
|---|
| 156 | integer, intent(in) :: nz, j |
|---|
| 157 | real, dimension(nz), intent(in) :: z, y |
|---|
| 158 | real, intent(out) :: dy_zj |
|---|
| 159 | |
|---|
| 160 | ! LOCAL VARIABLES |
|---|
| 161 | ! --------------- |
|---|
| 162 | real :: h1, h2, c1, c2, c3 |
|---|
| 163 | |
|---|
| 164 | ! CODE |
|---|
| 165 | ! ---- |
|---|
| 166 | if (j < 1 .or. j > nz - 2) then |
|---|
| 167 | dy_zj = -1. |
|---|
| 168 | else |
|---|
| 169 | h1 = z(j + 1) - z(j) |
|---|
| 170 | h2 = z(j + 2)- z(j + 1) |
|---|
| 171 | c1 = -(2*h1 + h2)/(h1*(h1 + h2)) |
|---|
| 172 | c2 = (h1 + h2)/(h1*h2) |
|---|
| 173 | c3 = -h1/(h2*(h1 + h2)) |
|---|
| 174 | dy_zj = c1*y(j) + c2*y(j + 1) + c3*y(j + 2) |
|---|
| 175 | endif |
|---|
| 176 | |
|---|
| 177 | END SUBROUTINE deriv1_onesided |
|---|
| 178 | !======================================================================= |
|---|
| 179 | |
|---|
| 180 | !======================================================================= |
|---|
| 181 | PURE SUBROUTINE colint(y,z,nz,i1,i2,integral) |
|---|
| 182 | !----------------------------------------------------------------------- |
|---|
| 183 | ! NAME |
|---|
| 184 | ! colint |
|---|
| 185 | ! |
|---|
| 186 | ! DESCRIPTION |
|---|
| 187 | ! Column integrates y on irregular grid. |
|---|
| 188 | ! |
|---|
| 189 | ! AUTHORS & DATE |
|---|
| 190 | ! N.S (Icarus 2010) |
|---|
| 191 | ! JB Clement, 2023-2025 |
|---|
| 192 | ! |
|---|
| 193 | ! NOTES |
|---|
| 194 | ! |
|---|
| 195 | !----------------------------------------------------------------------- |
|---|
| 196 | |
|---|
| 197 | ! DECLARATION |
|---|
| 198 | ! ----------- |
|---|
| 199 | implicit none |
|---|
| 200 | |
|---|
| 201 | ! ARGUMENTS |
|---|
| 202 | ! --------- |
|---|
| 203 | integer, intent(in) :: nz, i1, i2 |
|---|
| 204 | real, dimension(nz), intent(in) :: y, z |
|---|
| 205 | real, intent(out) :: integral |
|---|
| 206 | |
|---|
| 207 | ! LOCAL VARIABLES |
|---|
| 208 | ! --------------- |
|---|
| 209 | integer :: i |
|---|
| 210 | real, dimension(nz) :: dz |
|---|
| 211 | |
|---|
| 212 | ! CODE |
|---|
| 213 | ! ---- |
|---|
| 214 | dz(1) = (z(2) - 0.)/2 |
|---|
| 215 | do i = 2,nz - 1 |
|---|
| 216 | dz(i) = (z(i + 1) - z(i - 1))/2. |
|---|
| 217 | enddo |
|---|
| 218 | dz(nz) = z(nz) - z(nz - 1) |
|---|
| 219 | integral = sum(y(i1:i2)*dz(i1:i2)) |
|---|
| 220 | |
|---|
| 221 | END SUBROUTINE colint |
|---|
| 222 | !======================================================================= |
|---|
| 223 | |
|---|
| 224 | !======================================================================= |
|---|
| 225 | SUBROUTINE findroot(y1,y2,z1,z2,zr) |
|---|
| 226 | !----------------------------------------------------------------------- |
|---|
| 227 | ! NAME |
|---|
| 228 | ! findroot |
|---|
| 229 | ! |
|---|
| 230 | ! DESCRIPTION |
|---|
| 231 | ! Compute the root zr, between two values y1 and y2 at depth z1,z2. |
|---|
| 232 | ! |
|---|
| 233 | ! AUTHORS & DATE |
|---|
| 234 | ! L. Lange |
|---|
| 235 | ! JB Clement, 2023-2025 |
|---|
| 236 | ! |
|---|
| 237 | ! NOTES |
|---|
| 238 | ! |
|---|
| 239 | !----------------------------------------------------------------------- |
|---|
| 240 | |
|---|
| 241 | ! DECLARATION |
|---|
| 242 | ! ----------- |
|---|
| 243 | implicit none |
|---|
| 244 | |
|---|
| 245 | ! ARGUMENTS |
|---|
| 246 | ! --------- |
|---|
| 247 | real, intent(in) :: y1, y2 |
|---|
| 248 | real, intent(in) :: z1, z2 |
|---|
| 249 | real, intent(out) :: zr |
|---|
| 250 | |
|---|
| 251 | ! CODE |
|---|
| 252 | ! ---- |
|---|
| 253 | zr = (y1*z2 - y2*z1)/(y1 - y2) |
|---|
| 254 | |
|---|
| 255 | END SUBROUTINE findroot |
|---|
| 256 | !======================================================================= |
|---|
| 257 | |
|---|
| 258 | !======================================================================= |
|---|
| 259 | SUBROUTINE solve_tridiag(a,b,c,d,n,x,error) |
|---|
| 260 | !----------------------------------------------------------------------- |
|---|
| 261 | ! NAME |
|---|
| 262 | ! solve_tridiag |
|---|
| 263 | ! |
|---|
| 264 | ! DESCRIPTION |
|---|
| 265 | ! Solve a tridiagonal system Ax = d using the Thomas' algorithm |
|---|
| 266 | ! a: sub-diagonal |
|---|
| 267 | ! b: main diagonal |
|---|
| 268 | ! c: super-diagonal |
|---|
| 269 | ! d: right-hand side |
|---|
| 270 | ! x: solution |
|---|
| 271 | ! |
|---|
| 272 | ! AUTHORS & DATE |
|---|
| 273 | ! JB Clement, 2025 |
|---|
| 274 | ! |
|---|
| 275 | ! NOTES |
|---|
| 276 | ! |
|---|
| 277 | !----------------------------------------------------------------------- |
|---|
| 278 | |
|---|
| 279 | ! DECLARATION |
|---|
| 280 | ! ----------- |
|---|
| 281 | implicit none |
|---|
| 282 | |
|---|
| 283 | ! ARGUMENTS |
|---|
| 284 | ! --------- |
|---|
| 285 | integer, intent(in) :: n |
|---|
| 286 | real, dimension(n), intent(in) :: b, d |
|---|
| 287 | real, dimension(n - 1), intent(in) :: a, c |
|---|
| 288 | real, dimension(n), intent(out) :: x |
|---|
| 289 | integer, intent(out) :: error |
|---|
| 290 | |
|---|
| 291 | ! LOCAL VARIABLES |
|---|
| 292 | ! --------------- |
|---|
| 293 | integer :: i |
|---|
| 294 | real :: m |
|---|
| 295 | real, dimension(n) :: cp, dp |
|---|
| 296 | |
|---|
| 297 | ! CODE |
|---|
| 298 | ! ---- |
|---|
| 299 | ! Check stability: diagonally dominant condition |
|---|
| 300 | error = 0 |
|---|
| 301 | if (abs(b(1)) < abs(c(1))) then |
|---|
| 302 | error = 1 |
|---|
| 303 | return |
|---|
| 304 | endif |
|---|
| 305 | do i = 2,n - 1 |
|---|
| 306 | if (abs(b(i)) < abs(a(i - 1)) + abs(c(i))) then |
|---|
| 307 | error = 1 |
|---|
| 308 | return |
|---|
| 309 | endif |
|---|
| 310 | enddo |
|---|
| 311 | if (abs(b(n)) < abs(a(n - 1))) then |
|---|
| 312 | error = 1 |
|---|
| 313 | return |
|---|
| 314 | endif |
|---|
| 315 | |
|---|
| 316 | ! Initialization |
|---|
| 317 | cp(1) = c(1)/b(1) |
|---|
| 318 | dp(1) = d(1)/b(1) |
|---|
| 319 | |
|---|
| 320 | ! Forward sweep |
|---|
| 321 | do i = 2,n - 1 |
|---|
| 322 | m = b(i) - a(i - 1)*cp(i - 1) |
|---|
| 323 | cp(i) = c(i)/m |
|---|
| 324 | dp(i) = (d(i) - a(i - 1)*dp(i - 1))/m |
|---|
| 325 | enddo |
|---|
| 326 | m = b(n) - a(n - 1)*cp(n - 1) |
|---|
| 327 | dp(n) = (d(n) - a(n - 1)*dp(n - 1))/m |
|---|
| 328 | |
|---|
| 329 | ! Backward substitution |
|---|
| 330 | x(n) = dp(n) |
|---|
| 331 | do i = n - 1,1,-1 |
|---|
| 332 | x(i) = dp(i) - cp(i)*x(i + 1) |
|---|
| 333 | enddo |
|---|
| 334 | |
|---|
| 335 | END SUBROUTINE solve_tridiag |
|---|
| 336 | !======================================================================= |
|---|
| 337 | |
|---|
| 338 | !======================================================================= |
|---|
| 339 | SUBROUTINE solve_steady_heat(n,z,mz,kappa,mkappa,T_left,q_right,T) |
|---|
| 340 | !----------------------------------------------------------------------- |
|---|
| 341 | ! NAME |
|---|
| 342 | ! solve_steady_heat |
|---|
| 343 | ! |
|---|
| 344 | ! DESCRIPTION |
|---|
| 345 | ! Solve 1D steady-state heat equation with space-dependent thermal diffusivity. |
|---|
| 346 | ! Uses Thomas algorithm to solve tridiagonal system. |
|---|
| 347 | ! Left boundary: prescribed temperature (Dirichlet). |
|---|
| 348 | ! Right boundary: prescribed thermal flux (Neumann). |
|---|
| 349 | ! |
|---|
| 350 | ! AUTHORS & DATE |
|---|
| 351 | ! JB Clement, 2025 |
|---|
| 352 | ! |
|---|
| 353 | ! NOTES |
|---|
| 354 | ! Grid points at z, mid-grid points at mz. |
|---|
| 355 | ! Thermal diffusivity specified at both grids (kappa, mkappa). |
|---|
| 356 | !----------------------------------------------------------------------- |
|---|
| 357 | |
|---|
| 358 | ! DEPENDENCIES |
|---|
| 359 | ! ------------ |
|---|
| 360 | use stop_pem, only: stop_clean |
|---|
| 361 | |
|---|
| 362 | ! DECLARATION |
|---|
| 363 | ! ----------- |
|---|
| 364 | implicit none |
|---|
| 365 | |
|---|
| 366 | ! ARGUMENTS |
|---|
| 367 | ! --------- |
|---|
| 368 | integer, intent(in) :: n |
|---|
| 369 | real, dimension(n), intent(in) :: z, kappa |
|---|
| 370 | real, dimension(n - 1), intent(in) :: mz, mkappa |
|---|
| 371 | real, intent(in) :: T_left, q_right |
|---|
| 372 | real, dimension(n), intent(out) :: T |
|---|
| 373 | |
|---|
| 374 | ! LOCAL VARIABLES |
|---|
| 375 | ! --------------- |
|---|
| 376 | integer :: i, error |
|---|
| 377 | real, dimension(n) :: b, d |
|---|
| 378 | real, dimension(n - 1) :: a, c |
|---|
| 379 | |
|---|
| 380 | ! CODE |
|---|
| 381 | ! ---- |
|---|
| 382 | ! Initialization |
|---|
| 383 | a = 0.; b = 0.; c = 0.; d = 0. |
|---|
| 384 | |
|---|
| 385 | ! Left boundary condition (Dirichlet: prescribed temperature) |
|---|
| 386 | b(1) = 1. |
|---|
| 387 | d(1) = T_left |
|---|
| 388 | |
|---|
| 389 | ! Internal points |
|---|
| 390 | do i = 2,n - 1 |
|---|
| 391 | a(i - 1) = -mkappa(i - 1)/((mz(i) - mz(i - 1))*(z(i) - z(i - 1))) |
|---|
| 392 | c(i) = -mkappa(i)/((mz(i) - mz(i - 1))*(z(i + 1) - z(i))) |
|---|
| 393 | b(i) = -(a(i - 1) + c(i)) |
|---|
| 394 | enddo |
|---|
| 395 | |
|---|
| 396 | ! Right boundary condition (Neumann: prescribed temperature) |
|---|
| 397 | a(n - 1) = kappa(n - 1)/(z(n) - z(n - 1)) |
|---|
| 398 | b(n) = -kappa(n)/(z(n) - z(n - 1)) |
|---|
| 399 | d(n) = q_right |
|---|
| 400 | |
|---|
| 401 | ! Solve the tridiagonal linear system with the Thomas' algorithm |
|---|
| 402 | call solve_tridiag(a,b,c,d,n,T,error) |
|---|
| 403 | if (error /= 0) call stop_clean("solve_steady_heat","Unstable solving!",1) |
|---|
| 404 | |
|---|
| 405 | END SUBROUTINE solve_steady_heat |
|---|
| 406 | !======================================================================= |
|---|
| 407 | |
|---|
| 408 | end module math_toolkit |
|---|