p3m_functions module reference

Cubic interpolation functions.

More…

Functions/Subroutines

p3m_interpolation()

Set up a piecewise cubic interpolation from cell averages and estimated edge slopes and values.

p3m_limiter()

Adust a piecewise cubic reconstruction with a limiter that adjusts the edge values and slopes.

p3m_boundary_extrapolation()

Calculate the edge values and slopes at boundary cells as part of building a piecewise cubic sub-grid scale profiles.

build_cubic_interpolant()

Build cubic interpolant in cell k.

is_cubic_monotonic()

Check whether the cubic reconstruction in cell k is monotonic.

monotonize_cubic()

Monotonize a cubic curve by modifying the edge slopes.

Detailed Description

Date of creation: 2008.06.09 L. White.

This module contains p3m interpolation routines.

p3m interpolation is performed by estimating the edge values and slopes and constructing a cubic polynomial. We then make sure that the edge values are bounded and continuous and we then modify the slopes to get a monotonic cubic curve.

Function/Subroutine Documentation

subroutine p3m_functions/p3m_interpolation(N, h, u, edge_values, ppoly_S, ppoly_coef, h_neglect, answer_date)

Set up a piecewise cubic interpolation from cell averages and estimated edge slopes and values.

Cubic interpolation between edges.

The edge values and slopes MUST have been estimated prior to calling this routine.

It is assumed that the size of the array ‘u’ is equal to the number of cells defining ‘grid’ and ‘ppoly’. No consistency check is performed here.

Parameters:
  • n :: [in] Number of cells

  • h :: [in] cell widths (size N) [H]

  • u :: [in] cell averages (size N) in arbitrary units [A]

  • edge_values :: [inout] Edge value of polynomial [A]

  • ppoly_s :: [inout] Edge slope of polynomial [A H-1].

  • ppoly_coef :: [inout] Coefficients of polynomial [A]

  • h_neglect :: [in] A negligibly small width for the purpose of cell reconstructions [H]

  • answer_date :: [in] The vintage of the expressions to use

Call to:

p3m_limiter

Called from:

regrid_interp::regridding_set_ppolys

subroutine p3m_functions/p3m_limiter(N, h, u, edge_values, ppoly_S, ppoly_coef, h_neglect, answer_date)

Adust a piecewise cubic reconstruction with a limiter that adjusts the edge values and slopes.

The p3m limiter operates as follows:

  1. Edge values are bounded

  2. Discontinuous edge values are systematically averaged

  3. Loop on cells and do the following a. Build cubic curve b. Check if cubic curve is monotonic c. If not, monotonize cubic curve and rebuild it

Step 3 of the monotonization process leaves all edge values unchanged.

Parameters:
  • n :: [in] Number of cells

  • h :: [in] cell widths (size N) [H]

  • u :: [in] cell averages (size N) in arbitrary units [A]

  • edge_values :: [inout] Edge value of polynomial [A]

  • ppoly_s :: [inout] Edge slope of polynomial [A H-1]

  • ppoly_coef :: [inout] Coefficients of polynomial [A]

  • h_neglect :: [in] A negligibly small width for the purpose of cell reconstructions [H]

  • answer_date :: [in] The vintage of the expressions to use

Call to:

regrid_edge_values::average_discontinuous_edge_values regrid_edge_values::bound_edge_values build_cubic_interpolant hneglect_dflt is_cubic_monotonic monotonize_cubic

Called from:

p3m_interpolation

subroutine p3m_functions/p3m_boundary_extrapolation(N, h, u, edge_values, ppoly_S, ppoly_coef, h_neglect, h_neglect_edge)

Calculate the edge values and slopes at boundary cells as part of building a piecewise cubic sub-grid scale profiles.

The following explanations apply to the left boundary cell. The same reasoning holds for the right boundary cell.

A cubic needs to be built in the cell and requires four degrees of freedom, which are the edge values and slopes. The right edge values and slopes are taken to be that of the neighboring cell (i.e., the left edge value and slope of the neighboring cell). The left edge value and slope are determined by computing the parabola based on the cell average and the right edge value and slope. The resulting cubic is not necessarily monotonic and the slopes are subsequently modified to yield a monotonic cubic.

Parameters:
  • n :: [in] Number of cells

  • h :: [in] cell widths (size N) [H]

  • u :: [in] cell averages (size N) in arbitrary units [A]

  • edge_values :: [inout] Edge value of polynomial [A]

  • ppoly_s :: [inout] Edge slope of polynomial [A H-1]

  • ppoly_coef :: [inout] Coefficients of polynomial [A]

  • h_neglect :: [in] A negligibly small width for the purpose of cell reconstructions [H]

  • h_neglect_edge :: [in] A negligibly small width for the purpose of finding edge values [H]

Call to:

build_cubic_interpolant hneglect_dflt hneglect_edge_dflt is_cubic_monotonic monotonize_cubic

Called from:

regrid_interp::regridding_set_ppolys

subroutine p3m_functions/build_cubic_interpolant(h, k, edge_values, ppoly_S, ppoly_coef)

Build cubic interpolant in cell k.

Given edge values and edge slopes, compute coefficients of cubic in cell k.

NOTE: edge values and slopes MUST have been properly calculated prior to calling this routine.

Parameters:
  • h :: [in] cell widths (size N) [H]

  • k :: [in] The index of the cell to work on

  • edge_values :: [in] Edge value of polynomial in arbitrary units [A]

  • ppoly_s :: [in] Edge slope of polynomial [A H-1]

  • ppoly_coef :: [inout] Coefficients of polynomial [A]

Called from:

p3m_boundary_extrapolation p3m_limiter

function p3m_functions/is_cubic_monotonic(ppoly_coef, k) [logical]

Check whether the cubic reconstruction in cell k is monotonic.

This function checks whether the cubic curve in cell k is monotonic. If so, returns 1. Otherwise, returns 0.

The cubic is monotonic if the first derivative is single-signed in (0,1). Hence, we check whether the roots (if any) lie inside this interval. If there is no root or if both roots lie outside this interval, the cubic is monotonic.

Parameters:
  • ppoly_coef :: [in] Coefficients of cubic polynomial in arbitrary units [A]

  • k :: [in] The index of the cell to work on

Called from:

p3m_boundary_extrapolation p3m_limiter

subroutine p3m_functions/monotonize_cubic(h, u0_l, u0_r, sigma_l, sigma_r, slope, u1_l, u1_r)

Monotonize a cubic curve by modifying the edge slopes.

This routine takes care of monotonizing a cubic on [0,1] by modifying the edge slopes. The edge values are NOT modified. The cubic is entirely determined by the four degrees of freedom u0_l, u0_r, u1_l and u1_r.

u1_l and u1_r are the edge slopes expressed in the GLOBAL coordinate system.

The monotonization occurs as follows.

Parameters:
  • h :: [in] cell width [H]

  • u0_l :: [in] left edge value in arbitrary units [A]

  • u0_r :: [in] right edge value [A]

  • sigma_l :: [in] left 2nd-order slopes [A H-1]

  • sigma_r :: [in] right 2nd-order slopes [A H-1]

  • slope :: [in] limited PLM slope [A H-1]

  • u1_l :: [inout] left edge slopes [A H-1]

  • u1_r :: [inout] right edge slopes [A H-1]

Called from:

p3m_boundary_extrapolation p3m_limiter