MOM6
p3m_functions Module Reference

Functions/Subroutines

subroutine, public p3m_interpolation (N, h, u, ppoly_E, ppoly_S, ppoly_coefficients)
 
subroutine p3m_limiter (N, h, u, ppoly_E, ppoly_S, ppoly_coefficients)
 
subroutine, public p3m_boundary_extrapolation (N, h, u, ppoly_E, ppoly_S, ppoly_coefficients)
 
subroutine build_cubic_interpolant (h, k, ppoly_E, ppoly_S, ppoly_coefficients)
 
integer function is_cubic_monotonic (ppoly_coefficients, k)
 
subroutine monotonize_cubic (h, u0_l, u0_r, sigma_l, sigma_r, slope, u1_l, u1_r)
 

Variables

real, parameter h_neglect = 1.E-30
 

Function/Subroutine Documentation

◆ build_cubic_interpolant()

subroutine p3m_functions::build_cubic_interpolant ( real, dimension(:), intent(in)  h,
integer, intent(in)  k,
real, dimension(:,:), intent(in)  ppoly_E,
real, dimension(:,:), intent(in)  ppoly_S,
real, dimension(:,:), intent(inout)  ppoly_coefficients 
)
private

Definition at line 362 of file P3M_functions.F90.

Referenced by p3m_boundary_extrapolation(), and p3m_limiter().

362 !------------------------------------------------------------------------------
363 ! Given edge values and edge slopes, compute coefficients of cubic in cell k.
364 !
365 ! NOTE: edge values and slopes MUST have been properly calculated prior to
366 ! calling this routine.
367 !------------------------------------------------------------------------------
368 
369  ! Arguments
370  real, dimension(:), intent(in) :: h ! cell widths (size N)
371  integer, intent(in) :: k
372  real, dimension(:,:), intent(in) :: ppoly_e !Edge value of polynomial
373  real, dimension(:,:), intent(in) :: ppoly_s !Edge slope of polynomial
374  real, dimension(:,:), intent(inout) :: ppoly_coefficients !Coefficients of polynomial
375 
376  ! Local variables
377  real :: u0_l, u0_r ! edge values
378  real :: u1_l, u1_r ! edge slopes
379  real :: h_c ! cell width
380  real :: a0, a1, a2, a3 ! cubic coefficients
381 
382  h_c = h(k)
383 
384  u0_l = ppoly_e(k,1)
385  u0_r = ppoly_e(k,2)
386 
387  u1_l = ppoly_s(k,1) * h_c
388  u1_r = ppoly_s(k,2) * h_c
389 
390  a0 = u0_l
391  a1 = u1_l
392  a2 = 3.0 * ( u0_r - u0_l ) - u1_r - 2.0 * u1_l
393  a3 = u1_r + u1_l + 2.0 * ( u0_l - u0_r )
394 
395  ppoly_coefficients(k,1) = a0
396  ppoly_coefficients(k,2) = a1
397  ppoly_coefficients(k,3) = a2
398  ppoly_coefficients(k,4) = a3
399 
Here is the caller graph for this function:

◆ is_cubic_monotonic()

integer function p3m_functions::is_cubic_monotonic ( real, dimension(:,:), intent(in)  ppoly_coefficients,
integer, intent(in)  k 
)
private

Definition at line 407 of file P3M_functions.F90.

Referenced by p3m_boundary_extrapolation(), and p3m_limiter().

407 !------------------------------------------------------------------------------
408 ! This function checks whether the cubic curve in cell k is monotonic.
409 ! If so, returns 1. Otherwise, returns 0.
410 !
411 ! The cubic is monotonic if the first derivative is single-signed in [0,1].
412 ! Hence, we check whether the roots (if any) lie inside this interval. If there
413 ! is no root or if both roots lie outside this interval, the cubic is monotnic.
414 !------------------------------------------------------------------------------
415 
416  ! Arguments
417  real, dimension(:,:), intent(in) :: ppoly_coefficients
418  integer, intent(in) :: k
419 
420  ! Local variables
421  integer :: monotonic ! boolean indicating if monotonic or not
422  real :: a0, a1, a2, a3 ! cubic coefficients
423  real :: a, b, c ! coefficients of first derivative
424  real :: xi_0, xi_1 ! roots of first derivative (if any !)
425  real :: rho
426  real :: eps
427 
428  ! Define the radius of the ball around 0 and 1 in which all values are assumed
429  ! to be equal to 0 or 1, respectively
430  eps = 1e-14
431 
432  a0 = ppoly_coefficients(k,1)
433  a1 = ppoly_coefficients(k,2)
434  a2 = ppoly_coefficients(k,3)
435  a3 = ppoly_coefficients(k,4)
436 
437  a = a1
438  b = 2.0 * a2
439  c = 3.0 * a3
440 
441  xi_0 = -1.0
442  xi_1 = -1.0
443 
444  rho = b*b - 4.0*a*c
445 
446  if ( rho .GE. 0.0 ) then
447  if ( abs(c) .GT. 1e-15 ) then
448  xi_0 = 0.5 * ( -b - sqrt( rho ) ) / c
449  xi_1 = 0.5 * ( -b + sqrt( rho ) ) / c
450  else if ( abs(b) .GT. 1e-15 ) then
451  xi_0 = - a / b
452  xi_1 = - a / b
453  end if
454 
455  ! If one of the roots of the first derivative lies in (0,1),
456  ! the cubic is not monotonic.
457  if ( ( (xi_0 .GT. eps) .AND. (xi_0 .LT. 1.0-eps) ) .OR. &
458  ( (xi_1 .GT. eps) .AND. (xi_1 .LT. 1.0-eps) ) ) then
459  monotonic = 0
460  else
461  monotonic = 1
462  end if
463 
464  else ! there are no real roots --> cubic is monotonic
465  monotonic = 1
466  end if
467 
468  ! Set the return value
469  is_cubic_monotonic = monotonic
470 
Here is the caller graph for this function:

◆ monotonize_cubic()

subroutine p3m_functions::monotonize_cubic ( real, intent(in)  h,
real, intent(in)  u0_l,
real, intent(in)  u0_r,
real, intent(in)  sigma_l,
real, intent(in)  sigma_r,
real, intent(in)  slope,
real, intent(inout)  u1_l,
real, intent(inout)  u1_r 
)
private

Definition at line 478 of file P3M_functions.F90.

Referenced by p3m_boundary_extrapolation(), and p3m_limiter().

478 !------------------------------------------------------------------------------
479 ! This routine takes care of monotonizing a cubic on [0,1] by modifying the
480 ! edge slopes. The edge values are NOT modified. The cubic is entirely
481 ! determined by the four degrees of freedom u0_l, u0_r, u1_l and u1_r.
482 !
483 ! u1_l and u1_r are the edge slopes expressed in the GLOBAL coordinate system.
484 !
485 ! The monotonization occurs as follows.
486 
487 ! 1. The edge slopes are set to 0 if they are inconsistent with the limited
488 ! PLM slope
489 ! 2. We check whether we can find an inflexion point in [0,1]. At most one
490 ! inflexion point may exist.
491 ! (a) If there is no inflexion point, the cubic is monotonic.
492 ! (b) If there is one inflexion point and it lies outside [0,1], the
493 ! cubic is monotonic.
494 ! (c) If there is one inflexion point and it lies in [0,1] and the slope
495 ! at the location of the inflexion point is consistent, the cubic
496 ! is monotonic.
497 ! (d) If the inflexion point lies in [0,1] but the slope is inconsistent,
498 ! we go to (3) to shift the location of the inflexion point to the left
499 ! or to the right. To the left when the 2nd-order left slope is smaller
500 ! than the 2nd order right slope.
501 ! 3. Edge slopes are modified to shift the inflexion point, either onto the left
502 ! edge or onto the right edge.
503 !
504 !------------------------------------------------------------------------------
505 
506  ! Arguments
507  real, intent(in) :: h ! cell width
508  real, intent(in) :: u0_l, u0_r ! edge values
509  real, intent(in) :: sigma_l, sigma_r ! left and right 2nd-order slopes
510  real, intent(in) :: slope ! limited PLM slope
511  real, intent(inout) :: u1_l, u1_r ! edge slopes
512 
513  ! Local variables
514  integer :: found_ip
515  integer :: inflexion_l ! bool telling if inflex. pt must be on left
516  integer :: inflexion_r ! bool telling if inflex. pt must be on right
517  real :: eps
518  real :: a1, a2, a3
519  real :: u1_l_tmp ! trial left edge slope
520  real :: u1_r_tmp ! trial right edge slope
521  real :: xi_ip ! location of inflexion point
522  real :: slope_ip ! slope at inflexion point
523 
524  eps = 1e-14
525 
526  found_ip = 0
527  inflexion_l = 0
528  inflexion_r = 0
529 
530  ! If the edge slopes are inconsistent w.r.t. the limited PLM slope,
531  ! set them to zero
532  if ( u1_l*slope .LE. 0.0 ) then
533  u1_l = 0.0
534  end if
535 
536  if ( u1_r*slope .LE. 0.0 ) then
537  u1_r = 0.0
538  end if
539 
540  ! Compute the location of the inflexion point, which is the root
541  ! of the second derivative
542  a1 = h * u1_l
543  a2 = 3.0 * ( u0_r - u0_l ) - h*(u1_r + 2.0*u1_l)
544  a3 = h*(u1_r + u1_l) + 2.0*(u0_l - u0_r)
545 
546  ! There is a possible root (and inflexion point) only if a3 is nonzero.
547  ! When a3 is zero, the second derivative of the cubic is constant (the
548  ! cubic degenerates into a parabola) and no inflexion point exists.
549  if ( a3 .NE. 0.0 ) then
550  ! Location of inflexion point
551  xi_ip = - a2 / (3.0 * a3)
552 
553  ! If the inflexion point lies in [0,1], change boolean value
554  if ( (xi_ip .GE. 0.0) .AND. (xi_ip .LE. 1.0) ) then
555  found_ip = 1
556  end if
557  end if
558 
559  ! When there is an inflexion point within [0,1], check the slope
560  ! to see if it is consistent with the limited PLM slope. If not,
561  ! decide on which side we want to collapse the inflexion point.
562  ! If the inflexion point lies on one of the edges, the cubic is
563  ! guaranteed to be monotonic
564  if ( found_ip .EQ. 1 ) then
565  slope_ip = a1 + 2.0*a2*xi_ip + 3.0*a3*xi_ip*xi_ip
566 
567  ! Check whether slope is consistent
568  if ( slope_ip*slope .LT. 0.0 ) then
569  if ( abs(sigma_l) .LT. abs(sigma_r) ) then
570  inflexion_l = 1
571  else
572  inflexion_r = 1
573  end if
574  end if
575  end if ! found_ip
576 
577  ! At this point, if the cubic is not monotonic, we know where the
578  ! inflexion point should lie. When the cubic is monotonic, both
579  ! 'inflexion_l' and 'inflexion_r' are set to 0 and nothing is to be done.
580 
581  ! Move inflexion point on the left
582  if ( inflexion_l .EQ. 1 ) then
583 
584  u1_l_tmp = 1.5*(u0_r-u0_l)/h - 0.5*u1_r
585  u1_r_tmp = 3.0*(u0_r-u0_l)/h - 2.0*u1_l
586 
587  if ( (u1_l_tmp*slope .LT. 0.0) .AND. (u1_r_tmp*slope .LT. 0.0) ) then
588 
589  u1_l = 0.0
590  u1_r = 3.0 * (u0_r - u0_l) / h
591 
592  else if (u1_l_tmp*slope .LT. 0.0) then
593 
594  u1_r = u1_r_tmp
595  u1_l = 1.5*(u0_r - u0_l)/h - 0.5*u1_r
596 
597  else if (u1_r_tmp*slope .LT. 0.0) then
598 
599  u1_l = u1_l_tmp
600  u1_r = 3.0*(u0_r - u0_l)/h - 2.0*u1_l
601 
602  else
603 
604  u1_l = u1_l_tmp
605  u1_r = u1_r_tmp
606 
607  end if
608 
609  end if ! end treating case with inflexion point on the left
610 
611  ! Move inflexion point on the right
612  if ( inflexion_r .EQ. 1 ) then
613 
614  u1_l_tmp = 3.0*(u0_r-u0_l)/h - 2.0*u1_r
615  u1_r_tmp = 1.5*(u0_r-u0_l)/h - 0.5*u1_l
616 
617  if ( (u1_l_tmp*slope .LT. 0.0) .AND. (u1_r_tmp*slope .LT. 0.0) ) then
618 
619  u1_l = 3.0 * (u0_r - u0_l) / h
620  u1_r = 0.0
621 
622  else if (u1_l_tmp*slope .LT. 0.0) then
623 
624  u1_r = u1_r_tmp
625  u1_l = 3.0*(u0_r - u0_l)/h - 2.0*u1_r
626 
627  else if (u1_r_tmp*slope .LT. 0.0) then
628 
629  u1_l = u1_l_tmp
630  u1_r = 1.5*(u0_r - u0_l)/h - 0.5*u1_l
631 
632  else
633 
634  u1_l = u1_l_tmp
635  u1_r = u1_r_tmp
636 
637  end if
638 
639  end if ! end treating case with inflexion point on the right
640 
641  if ( abs(u1_l*h) .LT. eps ) then
642  u1_l = 0.0
643  end if
644 
645  if ( abs(u1_r*h) .LT. eps ) then
646  u1_r = 0.0
647  end if
648 
Here is the caller graph for this function:

◆ p3m_boundary_extrapolation()

subroutine, public p3m_functions::p3m_boundary_extrapolation ( integer, intent(in)  N,
real, dimension(:), intent(in)  h,
real, dimension(:), intent(in)  u,
real, dimension(:,:), intent(inout)  ppoly_E,
real, dimension(:,:), intent(inout)  ppoly_S,
real, dimension(:,:), intent(inout)  ppoly_coefficients 
)

Definition at line 204 of file P3M_functions.F90.

References build_cubic_interpolant(), h_neglect, is_cubic_monotonic(), and monotonize_cubic().

Referenced by regrid_interp::regridding_set_ppolys().

204 !------------------------------------------------------------------------------
205 ! The following explanations apply to the left boundary cell. The same
206 ! reasoning holds for the right boundary cell.
207 !
208 ! A cubic needs to be built in the cell and requires four degrees of freedom,
209 ! which are the edge values and slopes. The right edge values and slopes are
210 ! taken to be that of the neighboring cell (i.e., the left edge value and slope
211 ! of the neighboring cell). The left edge value and slope are determined by
212 ! computing the parabola based on the cell average and the right edge value
213 ! and slope. The resulting cubic is not necessarily monotonic and the slopes
214 ! are subsequently modified to yield a monotonic cubic.
215 !------------------------------------------------------------------------------
216 
217  ! Arguments
218  integer, intent(in) :: n ! Number of cells
219  real, dimension(:), intent(in) :: h ! cell widths (size N)
220  real, dimension(:), intent(in) :: u ! cell averages (size N)
221  real, dimension(:,:), intent(inout) :: ppoly_e !Edge value of polynomial
222  real, dimension(:,:), intent(inout) :: ppoly_s !Edge slope of polynomial
223  real, dimension(:,:), intent(inout) :: ppoly_coefficients !Coefficients of polynomial
224 
225  ! Local variables
226  integer :: i0, i1
227  integer :: monotonic
228  real :: u0, u1
229  real :: h0, h1
230  real :: b, c, d
231  real :: u0_l, u0_r
232  real :: u1_l, u1_r
233  real :: eps
234  real :: slope
235 
236  eps = 1e-10
237 
238  ! ----- Left boundary -----
239  i0 = 1
240  i1 = 2
241  h0 = h(i0) + eps
242  h1 = h(i1) + eps
243  u0 = u(i0)
244  u1 = u(i1)
245 
246  ! Compute the left edge slope in neighboring cell and express it in
247  ! the global coordinate system
248  b = ppoly_coefficients(i1,2)
249  u1_r = b / h1 ! derivative evaluated at xi = 0.0, expressed w.r.t. x
250 
251  ! Limit the right slope by the PLM limited slope
252  slope = 2.0 * ( u1 - u0 ) / ( h0 + h_neglect )
253  if ( abs(u1_r) .GT. abs(slope) ) then
254  u1_r = slope
255  end if
256 
257  ! The right edge value in the boundary cell is taken to be the left
258  ! edge value in the neighboring cell
259  u0_r = ppoly_e(i1,1)
260 
261  ! Given the right edge value and slope, we determine the left
262  ! edge value and slope by computing the parabola as determined by
263  ! the right edge value and slope and the boundary cell average
264  u0_l = 3.0 * u0 + 0.5 * h0*u1_r - 2.0 * u0_r
265  u1_l = ( - 6.0 * u0 - 2.0 * h0*u1_r + 6.0 * u0_r) / ( h0 + h_neglect )
266 
267  ! Check whether the edge values are monotonic. For example, if the left edge
268  ! value is larger than the right edge value while the slope is positive, the
269  ! edge values are inconsistent and we need to modify the left edge value
270  if ( (u0_r-u0_l) * slope .LT. 0.0 ) then
271  u0_l = u0_r
272  u1_l = 0.0
273  u1_r = 0.0
274  end if
275 
276  ! Store edge values and slope, build cubic and check monotonicity
277  ppoly_e(i0,1) = u0_l
278  ppoly_e(i0,2) = u0_r
279  ppoly_s(i0,1) = u1_l
280  ppoly_s(i0,2) = u1_r
281 
282  ! Store edge values and slope, build cubic and check monotonicity
283  call build_cubic_interpolant( h, i0, ppoly_e, ppoly_s, ppoly_coefficients )
284  monotonic = is_cubic_monotonic( ppoly_coefficients, i0 )
285 
286  if ( monotonic .EQ. 0 ) then
287  call monotonize_cubic( h0, u0_l, u0_r, 0.0, slope, slope, u1_l, u1_r )
288 
289  ! Rebuild cubic after monotonization
290  ppoly_s(i0,1) = u1_l
291  ppoly_s(i0,2) = u1_r
292  call build_cubic_interpolant( h, i0, ppoly_e, ppoly_s, ppoly_coefficients )
293 
294  end if
295 
296  ! ----- Right boundary -----
297  i0 = n-1
298  i1 = n
299  h0 = h(i0) + eps
300  h1 = h(i1) + eps
301  u0 = u(i0)
302  u1 = u(i1)
303 
304  ! Compute the right edge slope in neighboring cell and express it in
305  ! the global coordinate system
306  b = ppoly_coefficients(i0,2)
307  c = ppoly_coefficients(i0,3)
308  d = ppoly_coefficients(i0,4)
309  u1_l = (b + 2*c + 3*d) / ( h0 + h_neglect ) ! derivative evaluated at xi = 1.0
310 
311  ! Limit the left slope by the PLM limited slope
312  slope = 2.0 * ( u1 - u0 ) / ( h1 + h_neglect )
313  if ( abs(u1_l) .GT. abs(slope) ) then
314  u1_l = slope
315  end if
316 
317  ! The left edge value in the boundary cell is taken to be the right
318  ! edge value in the neighboring cell
319  u0_l = ppoly_e(i0,2)
320 
321  ! Given the left edge value and slope, we determine the right
322  ! edge value and slope by computing the parabola as determined by
323  ! the left edge value and slope and the boundary cell average
324  u0_r = 3.0 * u1 - 0.5 * h1*u1_l - 2.0 * u0_l
325  u1_r = ( 6.0 * u1 - 2.0 * h1*u1_l - 6.0 * u0_l) / ( h1 + h_neglect )
326 
327  ! Check whether the edge values are monotonic. For example, if the right edge
328  ! value is smaller than the left edge value while the slope is positive, the
329  ! edge values are inconsistent and we need to modify the right edge value
330  if ( (u0_r-u0_l) * slope .LT. 0.0 ) then
331  u0_r = u0_l
332  u1_l = 0.0
333  u1_r = 0.0
334  end if
335 
336  ! Store edge values and slope, build cubic and check monotonicity
337  ppoly_e(i1,1) = u0_l
338  ppoly_e(i1,2) = u0_r
339  ppoly_s(i1,1) = u1_l
340  ppoly_s(i1,2) = u1_r
341 
342  call build_cubic_interpolant( h, i1, ppoly_e, ppoly_s, ppoly_coefficients )
343  monotonic = is_cubic_monotonic( ppoly_coefficients, i1 )
344 
345  if ( monotonic .EQ. 0 ) then
346  call monotonize_cubic( h1, u0_l, u0_r, slope, 0.0, slope, u1_l, u1_r )
347 
348  ! Rebuild cubic after monotonization
349  ppoly_s(i1,1) = u1_l
350  ppoly_s(i1,2) = u1_r
351  call build_cubic_interpolant( h, i1, ppoly_e, ppoly_s, ppoly_coefficients )
352 
353  end if
354 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ p3m_interpolation()

subroutine, public p3m_functions::p3m_interpolation ( integer, intent(in)  N,
real, dimension(:), intent(in)  h,
real, dimension(:), intent(in)  u,
real, dimension(:,:), intent(inout)  ppoly_E,
real, dimension(:,:), intent(inout)  ppoly_S,
real, dimension(:,:), intent(inout)  ppoly_coefficients 
)

Definition at line 32 of file P3M_functions.F90.

References p3m_limiter().

Referenced by regrid_interp::regridding_set_ppolys().

32 !------------------------------------------------------------------------------
33 ! Cubic interpolation between edges.
34 !
35 ! The edge values and slopes MUST have been estimated prior to calling
36 ! this routine.
37 !
38 ! It is assumed that the size of the array 'u' is equal to the number of cells
39 ! defining 'grid' and 'ppoly'. No consistency check is performed here.
40 !------------------------------------------------------------------------------
41 
42  ! Arguments
43  integer, intent(in) :: n ! Number of cells
44  real, dimension(:), intent(in) :: h ! cell widths (size N)
45  real, dimension(:), intent(in) :: u ! cell averages (size N)
46  real, dimension(:,:), intent(inout) :: ppoly_e !Edge value of polynomial
47  real, dimension(:,:), intent(inout) :: ppoly_s !Edge slope of polynomial
48  real, dimension(:,:), intent(inout) :: ppoly_coefficients !Coefficients of polynomial
49 
50 
51  ! Call the limiter for p3m, which takes care of everything from
52  ! computing the coefficients of the cubic to monotonizing it.
53  ! This routine could be called directly instead of having to call
54  ! 'P3M_interpolation' first but we do that to provide an homogeneous
55  ! interface.
56 
57  call p3m_limiter( n, h, u, ppoly_e, ppoly_s, ppoly_coefficients )
58 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ p3m_limiter()

subroutine p3m_functions::p3m_limiter ( integer, intent(in)  N,
real, dimension(:), intent(in)  h,
real, dimension(:), intent(in)  u,
real, dimension(:,:), intent(inout)  ppoly_E,
real, dimension(:,:), intent(inout)  ppoly_S,
real, dimension(:,:), intent(inout)  ppoly_coefficients 
)
private

Definition at line 66 of file P3M_functions.F90.

References regrid_edge_values::average_discontinuous_edge_values(), regrid_edge_values::bound_edge_values(), build_cubic_interpolant(), h_neglect, is_cubic_monotonic(), and monotonize_cubic().

Referenced by p3m_interpolation().

66 !------------------------------------------------------------------------------
67 ! The p3m limiter operates as follows:
68 !
69 ! (1) Edge values are bounded
70 ! (2) Discontinuous edge values are systematically averaged
71 ! (3) Loop on cells and do the following
72 ! (a) Build cubic curve
73 ! (b) Check if cubic curve is monotonic
74 ! (c) If not, monotonize cubic curve and rebuild it
75 !
76 ! Step (3) of the monotonization process leaves all edge values unchanged.
77 !------------------------------------------------------------------------------
78 
79  ! Arguments
80  integer, intent(in) :: n ! Number of cells
81  real, dimension(:), intent(in) :: h ! cell widths (size N)
82  real, dimension(:), intent(in) :: u ! cell averages (size N)
83  real, dimension(:,:), intent(inout) :: ppoly_e !Edge value of polynomial
84  real, dimension(:,:), intent(inout) :: ppoly_s !Edge slope of polynomial
85  real, dimension(:,:), intent(inout) :: ppoly_coefficients !Coefficients of polynomial
86 
87 ! real, dimension(:,:), intent(inout) :: ppoly_coefficients
88 
89  ! Local variables
90  integer :: k ! loop index
91  integer :: monotonic ! boolean indicating whether the cubic is monotonic
92  real :: u0_l, u0_r ! edge values
93  real :: u1_l, u1_r ! edge slopes
94  real :: u_l, u_c, u_r ! left, center and right cell averages
95  real :: h_l, h_c, h_r ! left, center and right cell widths
96  real :: sigma_l, sigma_c, sigma_r ! left, center and right
97  ! van Leer slopes
98  real :: slope ! retained PLM slope
99  real :: eps
100 
101  eps = 1e-10
102 
103  ! 1. Bound edge values (boundary cells are assumed to be local extrema)
104  call bound_edge_values( n, h, u, ppoly_e )
105 
106  ! 2. Systematically average discontinuous edge values
107  call average_discontinuous_edge_values( n, ppoly_e )
108 
109 
110  ! 3. Loop on cells and do the following
111  ! (a) Build cubic curve
112  ! (b) Check if cubic curve is monotonic
113  ! (c) If not, monotonize cubic curve and rebuild it
114  do k = 1,n
115 
116  ! Get edge values, edge slopes and cell width
117  u0_l = ppoly_e(k,1)
118  u0_r = ppoly_e(k,2)
119  u1_l = ppoly_s(k,1)
120  u1_r = ppoly_s(k,2)
121 
122  ! Get cell widths and cell averages (boundary cells are assumed to
123  ! be local extrema for the sake of slopes)
124  u_c = u(k)
125  h_c = h(k)
126 
127  if ( k .EQ. 1 ) then
128  h_l = h(k)
129  u_l = u(k)
130  else
131  h_l = h(k-1)
132  u_l = u(k-1)
133  end if
134 
135  if ( k .EQ. n ) then
136  h_r = h(k)
137  u_r = u(k)
138  else
139  h_r = h(k+1)
140  u_r = u(k+1)
141  end if
142 
143  ! Compute limited slope
144  sigma_l = 2.0 * ( u_c - u_l ) / ( h_c + h_neglect )
145  sigma_c = 2.0 * ( u_r - u_l ) / ( h_l + 2.0*h_c + h_r + h_neglect )
146  sigma_r = 2.0 * ( u_r - u_c ) / ( h_c + h_neglect )
147 
148  if ( (sigma_l * sigma_r) .GT. 0.0 ) then
149  slope = sign( min(abs(sigma_l),abs(sigma_c),abs(sigma_r)), sigma_c )
150  else
151  slope = 0.0
152  end if
153 
154  ! If the slopes are close to zero in machine precision and in absolute
155  ! value, we set the slope to zero. This prevents asymmetric representation
156  ! near extrema.
157  if ( abs(u1_l*h_c) .LT. eps ) then
158  u1_l = 0.0
159  end if
160 
161  if ( abs(u1_r*h_c) .LT. eps ) then
162  u1_r = 0.0
163  end if
164 
165  ! The edge slopes are limited from above by the respective
166  ! one-sided slopes
167  if ( abs(u1_l) .GT. abs(sigma_l) ) then
168  u1_l = sigma_l
169  end if
170 
171  if ( abs(u1_r) .GT. abs(sigma_r) ) then
172  u1_r = sigma_r
173  end if
174 
175  ! Build cubic interpolant (compute the coefficients)
176  call build_cubic_interpolant( h, k, ppoly_e, ppoly_s, ppoly_coefficients )
177 
178  ! Check whether cubic is monotonic
179  monotonic = is_cubic_monotonic( ppoly_coefficients, k )
180 
181  ! If cubic is not monotonic, monotonize it by modifiying the
182  ! edge slopes, store the new edge slopes and recompute the
183  ! cubic coefficients
184  if ( monotonic .EQ. 0 ) then
185  call monotonize_cubic( h_c, u0_l, u0_r, sigma_l, sigma_r, slope, u1_l, u1_r )
186  end if
187 
188  ! Store edge slopes
189  ppoly_s(k,1) = u1_l
190  ppoly_s(k,2) = u1_r
191 
192  ! Recompute coefficients of cubic
193  call build_cubic_interpolant( h, k, ppoly_e, ppoly_s, ppoly_coefficients )
194 
195  end do ! loop on cells
196 
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ h_neglect

real, parameter p3m_functions::h_neglect = 1.E-30
private

Definition at line 24 of file P3M_functions.F90.

Referenced by p3m_boundary_extrapolation(), and p3m_limiter().

24 real, parameter :: h_neglect = 1.e-30