8 implicit none ;
private 10 public p3m_interpolation
11 public p3m_boundary_extrapolation
13 real,
parameter :: hneglect_dflt = 1.e-30
14 real,
parameter :: hneglect_edge_dflt = 1.e-10
28 subroutine p3m_interpolation( N, h, u, edge_values, ppoly_S, ppoly_coef, h_neglect, answers_2018 )
29 integer,
intent(in) :: N
30 real,
dimension(:),
intent(in) :: h
31 real,
dimension(:),
intent(in) :: u
32 real,
dimension(:,:),
intent(inout) :: edge_values
33 real,
dimension(:,:),
intent(inout) :: ppoly_S
34 real,
dimension(:,:),
intent(inout) :: ppoly_coef
35 real,
optional,
intent(in) :: h_neglect
37 logical,
optional,
intent(in) :: answers_2018
44 call p3m_limiter( n, h, u, edge_values, ppoly_s, ppoly_coef, h_neglect, answers_2018 )
46 end subroutine p3m_interpolation
61 subroutine p3m_limiter( N, h, u, edge_values, ppoly_S, ppoly_coef, h_neglect, answers_2018 )
62 integer,
intent(in) :: N
63 real,
dimension(:),
intent(in) :: h
64 real,
dimension(:),
intent(in) :: u
65 real,
dimension(:,:),
intent(inout) :: edge_values
66 real,
dimension(:,:),
intent(inout) :: ppoly_S
67 real,
dimension(:,:),
intent(inout) :: ppoly_coef
68 real,
optional,
intent(in) :: h_neglect
70 logical,
optional,
intent(in) :: answers_2018
79 real :: sigma_l, sigma_c, sigma_r
84 hneglect = hneglect_dflt ;
if (
present(h_neglect)) hneglect = h_neglect
89 call bound_edge_values( n, h, u, edge_values, hneglect, answers_2018 )
92 call average_discontinuous_edge_values( n, edge_values )
102 u0_l = edge_values(k,1)
103 u0_r = edge_values(k,2)
129 sigma_l = 2.0 * ( u_c - u_l ) / ( h_c + hneglect )
130 sigma_c = 2.0 * ( u_r - u_l ) / ( h_l + 2.0*h_c + h_r + hneglect )
131 sigma_r = 2.0 * ( u_r - u_c ) / ( h_c + hneglect )
133 if ( (sigma_l * sigma_r) > 0.0 )
then 134 slope = sign( min(abs(sigma_l),abs(sigma_c),abs(sigma_r)), sigma_c )
140 if ( abs(u1_l*h_c) < epsilon(u_c)*abs(u_c) ) u1_l = 0.0
141 if ( abs(u1_r*h_c) < epsilon(u_c)*abs(u_c) ) u1_r = 0.0
145 if ( abs(u1_l) > abs(sigma_l) )
then 149 if ( abs(u1_r) > abs(sigma_r) )
then 154 call build_cubic_interpolant( h, k, edge_values, ppoly_s, ppoly_coef )
157 monotonic = is_cubic_monotonic( ppoly_coef, k )
162 if ( .not.monotonic )
then 163 call monotonize_cubic( h_c, u0_l, u0_r, sigma_l, sigma_r, slope, u1_l, u1_r )
171 call build_cubic_interpolant( h, k, edge_values, ppoly_s, ppoly_coef )
175 end subroutine p3m_limiter
191 subroutine p3m_boundary_extrapolation( N, h, u, edge_values, ppoly_S, ppoly_coef, &
192 h_neglect, h_neglect_edge )
193 integer,
intent(in) :: N
194 real,
dimension(:),
intent(in) :: h
195 real,
dimension(:),
intent(in) :: u
196 real,
dimension(:,:),
intent(inout) :: edge_values
197 real,
dimension(:,:),
intent(inout) :: ppoly_S
198 real,
dimension(:,:),
intent(inout) :: ppoly_coef
199 real,
optional,
intent(in) :: h_neglect
201 real,
optional,
intent(in) :: h_neglect_edge
212 real :: hNeglect, hNeglect_edge
214 hneglect = hneglect_dflt ;
if (
present(h_neglect)) hneglect = h_neglect
215 hneglect_edge = hneglect_edge_dflt ;
if (
present(h_neglect_edge)) hneglect_edge = h_neglect_edge
220 h0 = h(i0) + hneglect_edge
221 h1 = h(i1) + hneglect_edge
231 slope = 2.0 * ( u1 - u0 ) / ( h0 + hneglect )
232 if ( abs(u1_r) > abs(slope) )
then 238 u0_r = edge_values(i1,1)
243 u0_l = 3.0 * u0 + 0.5 * h0*u1_r - 2.0 * u0_r
244 u1_l = ( - 6.0 * u0 - 2.0 * h0*u1_r + 6.0 * u0_r) / ( h0 + hneglect )
249 if ( (u0_r-u0_l) * slope < 0.0 )
then 256 edge_values(i0,1) = u0_l
257 edge_values(i0,2) = u0_r
262 call build_cubic_interpolant( h, i0, edge_values, ppoly_s, ppoly_coef )
263 monotonic = is_cubic_monotonic( ppoly_coef, i0 )
265 if ( .not.monotonic )
then 266 call monotonize_cubic( h0, u0_l, u0_r, 0.0, slope, slope, u1_l, u1_r )
271 call build_cubic_interpolant( h, i0, edge_values, ppoly_s, ppoly_coef )
278 h0 = h(i0) + hneglect_edge
279 h1 = h(i1) + hneglect_edge
288 u1_l = (b + 2*c + 3*d) / ( h0 + hneglect )
291 slope = 2.0 * ( u1 - u0 ) / ( h1 + hneglect )
292 if ( abs(u1_l) > abs(slope) )
then 298 u0_l = edge_values(i0,2)
303 u0_r = 3.0 * u1 - 0.5 * h1*u1_l - 2.0 * u0_l
304 u1_r = ( 6.0 * u1 - 2.0 * h1*u1_l - 6.0 * u0_l) / ( h1 + hneglect )
309 if ( (u0_r-u0_l) * slope < 0.0 )
then 316 edge_values(i1,1) = u0_l
317 edge_values(i1,2) = u0_r
321 call build_cubic_interpolant( h, i1, edge_values, ppoly_s, ppoly_coef )
322 monotonic = is_cubic_monotonic( ppoly_coef, i1 )
324 if ( .not.monotonic )
then 325 call monotonize_cubic( h1, u0_l, u0_r, slope, 0.0, slope, u1_l, u1_r )
330 call build_cubic_interpolant( h, i1, edge_values, ppoly_s, ppoly_coef )
334 end subroutine p3m_boundary_extrapolation
343 subroutine build_cubic_interpolant( h, k, edge_values, ppoly_S, ppoly_coef )
344 real,
dimension(:),
intent(in) :: h
345 integer,
intent(in) :: k
346 real,
dimension(:,:),
intent(in) :: edge_values
347 real,
dimension(:,:),
intent(in) :: ppoly_S
348 real,
dimension(:,:),
intent(inout) :: ppoly_coef
354 real :: a0, a1, a2, a3
358 u0_l = edge_values(k,1)
359 u0_r = edge_values(k,2)
361 u1_l = ppoly_s(k,1) * h_c
362 u1_r = ppoly_s(k,2) * h_c
366 a2 = 3.0 * ( u0_r - u0_l ) - u1_r - 2.0 * u1_l
367 a3 = u1_r + u1_l + 2.0 * ( u0_l - u0_r )
374 end subroutine build_cubic_interpolant
385 logical function is_cubic_monotonic( ppoly_coef, k )
386 real,
dimension(:,:),
intent(in) :: ppoly_coef
387 integer,
intent(in) :: k
392 b = 2.0 * ppoly_coef(k,3)
393 c = 3.0 * ppoly_coef(k,4)
396 if (b*b - 4.0*a*c <= 0.0)
then 397 is_cubic_monotonic = .true.
398 elseif (a * (a + (b + c)) < 0.0)
then 399 is_cubic_monotonic = .false.
400 elseif (b * (b + 2.0*c) < 0.0)
then 401 is_cubic_monotonic = .false.
403 is_cubic_monotonic = .true.
406 end function is_cubic_monotonic
435 subroutine monotonize_cubic( h, u0_l, u0_r, sigma_l, sigma_r, slope, u1_l, u1_r )
436 real,
intent(in) :: h
437 real,
intent(in) :: u0_l
438 real,
intent(in) :: u0_r
439 real,
intent(in) :: sigma_l
440 real,
intent(in) :: sigma_r
441 real,
intent(in) :: slope
442 real,
intent(inout) :: u1_l
443 real,
intent(inout) :: u1_r
446 logical :: inflexion_l
447 logical :: inflexion_r
455 inflexion_l = .false.
456 inflexion_r = .false.
460 if ( u1_l*slope <= 0.0 )
then 464 if ( u1_r*slope <= 0.0 )
then 471 a2 = 3.0 * ( u0_r - u0_l ) - h*(u1_r + 2.0*u1_l)
472 a3 = h*(u1_r + u1_l) + 2.0*(u0_l - u0_r)
477 if ( a3 /= 0.0 )
then 479 xi_ip = - a2 / (3.0 * a3)
482 if ( (xi_ip >= 0.0) .AND. (xi_ip <= 1.0) )
then 493 slope_ip = a1 + 2.0*a2*xi_ip + 3.0*a3*xi_ip*xi_ip
496 if ( slope_ip*slope < 0.0 )
then 497 if ( abs(sigma_l) < abs(sigma_r) )
then 510 if ( inflexion_l )
then 512 u1_l_tmp = 1.5*(u0_r-u0_l)/h - 0.5*u1_r
513 u1_r_tmp = 3.0*(u0_r-u0_l)/h - 2.0*u1_l
515 if ( (u1_l_tmp*slope < 0.0) .AND. (u1_r_tmp*slope < 0.0) )
then 518 u1_r = 3.0 * (u0_r - u0_l) / h
520 elseif (u1_l_tmp*slope < 0.0)
then 523 u1_l = 1.5*(u0_r - u0_l)/h - 0.5*u1_r
525 elseif (u1_r_tmp*slope < 0.0)
then 528 u1_r = 3.0*(u0_r - u0_l)/h - 2.0*u1_l
540 if ( inflexion_r )
then 542 u1_l_tmp = 3.0*(u0_r-u0_l)/h - 2.0*u1_r
543 u1_r_tmp = 1.5*(u0_r-u0_l)/h - 0.5*u1_l
545 if ( (u1_l_tmp*slope < 0.0) .AND. (u1_r_tmp*slope < 0.0) )
then 547 u1_l = 3.0 * (u0_r - u0_l) / h
550 elseif (u1_l_tmp*slope < 0.0)
then 553 u1_l = 3.0*(u0_r - u0_l)/h - 2.0*u1_r
555 elseif (u1_r_tmp*slope < 0.0)
then 558 u1_r = 1.5*(u0_r - u0_l)/h - 0.5*u1_l
570 if ( abs(u1_l*h) < epsilon(u0_l) * (abs(u0_l) + abs(u0_r)) ) u1_l = 0.0
571 if ( abs(u1_r*h) < epsilon(u0_l) * (abs(u0_l) + abs(u0_r)) ) u1_r = 0.0
573 end subroutine monotonize_cubic
Cubic interpolation functions.
Edge value estimation for high-order resconstruction.