plm_functions module reference

Piecewise linear reconstruction functions.

More…

Functions/Subroutines

plm_slope_wa()

Returns a limited PLM slope following White and Adcroft, 2008, in the same arbitrary units [A] as the input values.

plm_slope_cw()

Returns a limited PLM slope following Colella and Woodward 1984, in the same arbitrary units as the input values [A].

plm_monotonized_slope()

Returns a limited PLM slope following Colella and Woodward 1984, in the same arbitrary units as the input values [A].

plm_extrapolate_slope()

Returns a PLM slope using h2 extrapolation from a cell to the left, in the same arbitrary units as the input values [A].

plm_reconstruction()

Reconstruction by linear polynomials within each cell.

plm_boundary_extrapolation()

Reconstruction by linear polynomials within boundary cells.

Detailed Description

Date of creation: 2008.06.06 L. White.

This module contains routines that handle one-dimensional finite volume reconstruction using the piecewise linear method (PLM).

Function/Subroutine Documentation

function plm_functions/plm_slope_wa(h_l, h_c, h_r, h_neglect, u_l, u_c, u_r) [pure]

Returns a limited PLM slope following White and Adcroft, 2008, in the same arbitrary units [A] as the input values. Note that this is not the same as the Colella and Woodward method.

Parameters:
  • h_l :: [in] Thickness of left cell in arbitrary grid thickness units [H]

  • h_c :: [in] Thickness of center cell in arbitrary grid thickness units [H]

  • h_r :: [in] Thickness of right cell in arbitrary grid thickness units [H]

  • h_neglect :: [in] A negligible thickness [H]

  • u_l :: [in] Value of left cell in arbitrary units [A]

  • u_c :: [in] Value of center cell in arbitrary units [A]

  • u_r :: [in] Value of right cell in arbitrary units [A]

Called from:

mom_ale::ale_plm_edge_values plm_reconstruction

function plm_functions/plm_slope_cw(h_l, h_c, h_r, h_neglect, u_l, u_c, u_r) [pure]

Returns a limited PLM slope following Colella and Woodward 1984, in the same arbitrary units as the input values [A].

Parameters:
  • h_l :: [in] Thickness of left cell in arbitrary grid thickness units [H]

  • h_c :: [in] Thickness of center cell in arbitrary grid thickness units [H]

  • h_r :: [in] Thickness of right cell in arbitrary grid thickness units [H]

  • h_neglect :: [in] A negligible thickness [H]

  • u_l :: [in] Value of left cell in arbitrary units [A]

  • u_c :: [in] Value of center cell in arbitrary units [A]

  • u_r :: [in] Value of right cell in arbitrary units [A]

function plm_functions/plm_monotonized_slope(u_l, u_c, u_r, s_l, s_c, s_r) [pure]

Returns a limited PLM slope following Colella and Woodward 1984, in the same arbitrary units as the input values [A].

Parameters:
  • u_l :: [in] Value of left cell in arbitrary units [A]

  • u_c :: [in] Value of center cell in arbitrary units [A]

  • u_r :: [in] Value of right cell in arbitrary units [A]

  • s_l :: [in] PLM slope of left cell [A]

  • s_c :: [in] PLM slope of center cell [A]

  • s_r :: [in] PLM slope of right cell [A]

Called from:

mom_ale::ale_plm_edge_values plm_reconstruction

function plm_functions/plm_extrapolate_slope(h_l, h_c, h_neglect, u_l, u_c) [pure]

Returns a PLM slope using h2 extrapolation from a cell to the left, in the same arbitrary units as the input values [A]. Use the negative to extrapolate from the cell to the right.

Parameters:
  • h_l :: [in] Thickness of left cell in arbitrary grid thickness units [H]

  • h_c :: [in] Thickness of center cell in arbitrary grid thickness units [H]

  • h_neglect :: [in] A negligible thickness [H]

  • u_l :: [in] Value of left cell in arbitrary units [A]

  • u_c :: [in] Value of center cell in arbitrary units [A]

Called from:

mom_ale::ale_plm_edge_values plm_boundary_extrapolation

subroutine plm_functions/plm_reconstruction(N, h, u, edge_values, ppoly_coef, h_neglect)

Reconstruction by linear polynomials within each cell.

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 values of piecewise polynomials, with the same units as u [A].

  • ppoly_coef :: [inout] coefficients of piecewise polynomials, mainly with the same units as u [A].

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

Call to:

hneglect_dflt plm_monotonized_slope plm_slope_wa

Called from:

mom_remapping::build_reconstructions_1d regrid_interp::regridding_set_ppolys mom_remapping::remapping_unit_tests

subroutine plm_functions/plm_boundary_extrapolation(N, h, u, edge_values, ppoly_coef, h_neglect)

Reconstruction by linear polynomials within boundary cells.

The left and right edge values in the left and right boundary cells, respectively, are estimated using a linear extrapolation within the cells.

This extrapolation is EXACT when the underlying profile is linear.

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 values of piecewise polynomials, with the same units as u [A].

  • ppoly_coef :: [inout] coefficients of piecewise polynomials, mainly with the same units as u [A].

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

Call to:

hneglect_dflt plm_extrapolate_slope

Called from:

mom_remapping::build_reconstructions_1d regrid_interp::regridding_set_ppolys