PETSc version 3.15.5
Fix/Edit manual page

PetscLimiterLimit

Limit the flux

Synopsis

#include "petscfv.h" 
PetscErrorCode PetscLimiterLimit(PetscLimiter lim, PetscReal flim, PetscReal *phi)

Input Parameters

lim - The PetscLimiter
flim - The input field

Output Parameter

phi - The limited field

Note: Limiters given in symmetric form following Berger, Aftosmis, and Murman 2005

The classical flux-limited formulation is psi(r) where

r = (u[0] - u[-1]) / (u[1] - u[0])

The second order TVD region is bounded by

psi_minmod(r) = min(r,1)      and        psi_superbee(r) = min(2, 2r, max(1,r))

where all limiters are implicitly clipped to be non-negative. A more convenient slope-limited form is psi(r) =
phi(r)(r+1)/2 in which the reconstructed interface values are

u(v) = u[0] + phi(r) (grad u)[0] v

where v is the vector from centroid to quadrature point. In these variables, the usual limiters become

phi_minmod(r) = 2 min(1/(1+r),r/(1+r))   phi_superbee(r) = 2 min(2/(1+r), 2r/(1+r), max(1,r)/(1+r))

For a nicer symmetric formulation, rewrite in terms of

f = (u[0] - u[-1]) / (u[1] - u[-1])

where r(f) = f/(1-f). Not that r(1-f) = (1-f)/f = 1/r(f) so the symmetry condition

phi(r) = phi(1/r)

becomes

w(f) = w(1-f).

The limiters below implement this final form w(f). The reference methods are

w_minmod(f) = 2 min(f,(1-f))             w_superbee(r) = 4 min((1-f), f)

See Also

PetscLimiterSetType(), PetscLimiterCreate()

Level

beginner

Location

src/dm/dt/fv/interface/fv.c
Index of all FV routines
Table of Contents for all manual pages
Index of all manual pages