recon1d_ppm_h4_2019 module reference

Piecewise Parabolic Method 1D reconstruction with h4 interpolation for edges.

More…

Data Types

ppm_h4_2019

PPM reconstruction following White and Adcroft, 2008.

Functions/Subroutines

init()

Initialize a 1D PPM_H4_2019 reconstruction for n cells.

reconstruct()

Calculate a 1D PPM_H4_2019 reconstructions based on h(:) and u(:)

end_value_h4()

Determine a one-sided 4th order polynomial fit of u to the data points for the purposes of specifying edge values, as described in the appendix of White and Adcroft JCP 2008.

f()

Value of PPM_H4_2019 reconstruction at a point in cell k [A].

dfdx()

Derivative of PPM_H4_2019 reconstruction at a point in cell k [A].

average()

Average between xa and xb for cell k of a 1D PPM reconstruction [A].

destroy()

Deallocate the PPM_H4_2019 reconstruction.

check_reconstruction()

Checks the PPM_H4_2019 reconstruction for consistency.

unit_tests()

Runs PPM_H4_2019 reconstruction unit tests and returns True for any fails, False otherwise.

Detailed Description

This implementation of PPM follows White and Adcroft 2008 [87], with cells resorting to PCM for extrema including first and last cells in column. This scheme differs from Colella and Woodward, 1984 [17], in the method of first estimating the fourth-order accurate edge values. This uses numerical expressions refactored at the beginning of 2019. The first and last cells are always limited to PCM.

Type Documentation

type  recon1d_ppm_h4_2019/ppm_h4_2019

PPM reconstruction following White and Adcroft, 2008.

Type fields:
  • % ul :: real, dimension(:), allocatable, private Left edge value [A].

  • % ur :: real, dimension(:), allocatable, private Right edge value [A].

  • % init :: procedure, private Implementation of the PPM_H4_2019 initialization.

  • % reconstruct :: procedure, private Implementation of the PPM_H4_2019 reconstruction.

  • % average :: procedure, private Implementation of the PPM_H4_2019 average over an interval [A].

  • % f :: procedure, private Implementation of evaluating the PPM_H4_2019 reconstruction at a point [A].

  • % dfdx :: procedure, private Implementation of the derivative of the PPM_H4_2019 reconstruction at a point [A].

  • % destroy :: procedure, private Implementation of deallocation for PPM_H4_2019.

  • % check_reconstruction :: procedure, private Implementation of check reconstruction for the PPM_H4_2019 reconstruction.

  • % unit_tests :: procedure, private Implementation of unit tests for the PPM_H4_2019 reconstruction.

  • % init_parent :: procedure, private Duplicate interface to

  • % reconstruct_parent :: procedure, private Duplicate interface to

[source]

Function/Subroutine Documentation

subroutine recon1d_ppm_h4_2019/init(this, n, h_neglect, check)

Initialize a 1D PPM_H4_2019 reconstruction for n cells.

Parameters:
  • this :: this [out] This reconstruction

  • n :: n [in] Number of cells in this column

  • h_neglect :: h_neglect [in] A negligibly small width used in cell reconstructions [H]

  • check :: check [in] If true, enable some consistency checking

[source]

subroutine recon1d_ppm_h4_2019/reconstruct(this, h, u)

Calculate a 1D PPM_H4_2019 reconstructions based on h(:) and u(:)

Parameters:
  • this :: this [inout] This reconstruction

  • h :: h [in] Grid spacing (thickness) [typically H]

  • u :: u [in] Cell mean values [A]

Call to:

end_value_h4

[source]

subroutine recon1d_ppm_h4_2019/end_value_h4(dz, u, Csys)

Determine a one-sided 4th order polynomial fit of u to the data points for the purposes of specifying edge values, as described in the appendix of White and Adcroft JCP 2008.

Parameters:
  • dz :: dz [in] The thicknesses of 4 layers, starting at the edge [H]. The values of dz must be positive.

  • u :: u [in] The average properties of 4 layers, starting at the edge [A]

  • csys :: [out] The four coefficients of a 4th order polynomial fit of u as a function of z [A H-(n-1)]

Called from:

reconstruct

[source]

function  recon1d_ppm_h4_2019/f(this, k, x)

Value of PPM_H4_2019 reconstruction at a point in cell k [A].

Parameters:
  • this :: this [in] This reconstruction

  • k :: k [in] Cell number

  • x :: x [in] Non-dimensional position within element [nondim]

[source]

function  recon1d_ppm_h4_2019/dfdx(this, k, x)

Derivative of PPM_H4_2019 reconstruction at a point in cell k [A].

Parameters:
  • this :: this [in] This reconstruction

  • k :: k [in] Cell number

  • x :: x [in] Non-dimensional position within element [nondim]

[source]

function  recon1d_ppm_h4_2019/average(this, k, xa, xb)

Average between xa and xb for cell k of a 1D PPM reconstruction [A].

Parameters:
  • this :: this [in] This reconstruction

  • k :: k [in] Cell number

  • xa :: xa [in] Start of averaging interval on element (0 to 1)

  • xb :: xb [in] End of averaging interval on element (0 to 1)

[source]

subroutine recon1d_ppm_h4_2019/destroy(this)

Deallocate the PPM_H4_2019 reconstruction.

Parameters:

this :: this [inout] This reconstruction

[source]

function  recon1d_ppm_h4_2019/check_reconstruction(this, h, u)

Checks the PPM_H4_2019 reconstruction for consistency.

Parameters:
  • this :: this [in] This reconstruction

  • h :: h [in] Grid spacing (thickness) [typically H]

  • u :: u [in] Cell mean values [A]

[source]

function  recon1d_ppm_h4_2019/unit_tests(this, verbose, stdout, stderr)

Runs PPM_H4_2019 reconstruction unit tests and returns True for any fails, False otherwise.

Parameters:
  • this :: this [inout] This reconstruction

  • verbose :: verbose [in] True, if verbose

  • stdout :: stdout [in] I/O channel for stdout

  • stderr :: stderr [in] I/O channel for stderr

[source]

[source]