MFEM  v4.5.2
Finite element discretization library
Public Member Functions | Protected Attributes | List of all members
mfem::DeltaCoefficient Class Reference

Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent C-function. More...

#include <coefficient.hpp>

Inheritance diagram for mfem::DeltaCoefficient:
[legend]
Collaboration diagram for mfem::DeltaCoefficient:
[legend]

Public Member Functions

 DeltaCoefficient ()
 Construct a unit delta function centered at (0.0,0.0,0.0) More...
 
 DeltaCoefficient (double x, double s)
 Construct a delta function scaled by s and centered at (x,0.0,0.0) More...
 
 DeltaCoefficient (double x, double y, double s)
 Construct a delta function scaled by s and centered at (x,y,0.0) More...
 
 DeltaCoefficient (double x, double y, double z, double s)
 Construct a delta function scaled by s and centered at (x,y,z) More...
 
void SetTime (double t)
 Set the time for internally stored coefficients. More...
 
void SetDeltaCenter (const Vector &center)
 Set the center location of the delta function. More...
 
void SetScale (double s_)
 Set the scale value multiplying the delta function. More...
 
void SetFunction (double(*f)(double))
 Set a time-dependent function that multiplies the Scale(). More...
 
void SetTol (double tol_)
 Set the tolerance used during projection onto GridFunction to identify the Mesh vertex where the Center() of the delta function lies. (default 1e-12) More...
 
void SetWeight (Coefficient *w)
 Set a weight Coefficient that multiplies the DeltaCoefficient. More...
 
const double * Center ()
 
double Scale ()
 Return the scale factor times the optional time dependent function. Returns \( s T(t) \) with \( T(t) = 1 \) when not set by the user. More...
 
double Tol ()
 Return the tolerance used to identify the mesh vertices. More...
 
CoefficientWeight ()
 See SetWeight() for description of the weight Coefficient. More...
 
void GetDeltaCenter (Vector &center)
 Write the center of the delta function into center. More...
 
virtual double EvalDelta (ElementTransformation &T, const IntegrationPoint &ip)
 The value of the function assuming we are evaluating at the delta center. More...
 
virtual double Eval (ElementTransformation &T, const IntegrationPoint &ip)
 A DeltaFunction cannot be evaluated. Calling this method will cause an MFEM error, terminating the application. More...
 
virtual ~DeltaCoefficient ()
 
- Public Member Functions inherited from mfem::Coefficient
 Coefficient ()
 
double GetTime ()
 Get the time for time dependent coefficients. More...
 
double Eval (ElementTransformation &T, const IntegrationPoint &ip, double t)
 Evaluate the coefficient in the element described by T at the point ip at time t. More...
 
virtual void Project (QuadratureFunction &qf)
 Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points. More...
 
virtual ~Coefficient ()
 

Protected Attributes

double center [3]
 
double scale
 
double tol
 
Coefficientweight
 
int sdim
 
double(* tdf )(double)
 
- Protected Attributes inherited from mfem::Coefficient
double time
 

Detailed Description

Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent C-function.

\( F(x,t) = w(x,t) s T(t) d(x - xc) \)

where w is the optional weight coefficient, s is a scale factor T is an optional time-dependent function and d is a delta function.

WARNING this cannot be used as a normal coefficient. The usual Eval method is disabled.

Definition at line 334 of file coefficient.hpp.

Constructor & Destructor Documentation

◆ DeltaCoefficient() [1/4]

mfem::DeltaCoefficient::DeltaCoefficient ( )
inline

Construct a unit delta function centered at (0.0,0.0,0.0)

Definition at line 345 of file coefficient.hpp.

◆ DeltaCoefficient() [2/4]

mfem::DeltaCoefficient::DeltaCoefficient ( double  x,
double  s 
)
inline

Construct a delta function scaled by s and centered at (x,0.0,0.0)

Definition at line 352 of file coefficient.hpp.

◆ DeltaCoefficient() [3/4]

mfem::DeltaCoefficient::DeltaCoefficient ( double  x,
double  y,
double  s 
)
inline

Construct a delta function scaled by s and centered at (x,y,0.0)

Definition at line 359 of file coefficient.hpp.

◆ DeltaCoefficient() [4/4]

mfem::DeltaCoefficient::DeltaCoefficient ( double  x,
double  y,
double  z,
double  s 
)
inline

Construct a delta function scaled by s and centered at (x,y,z)

Definition at line 366 of file coefficient.hpp.

◆ ~DeltaCoefficient()

virtual mfem::DeltaCoefficient::~DeltaCoefficient ( )
inlinevirtual

Definition at line 421 of file coefficient.hpp.

Member Function Documentation

◆ Center()

const double* mfem::DeltaCoefficient::Center ( )
inline

Return a pointer to a c-array representing the center of the delta function.

Definition at line 399 of file coefficient.hpp.

◆ Eval()

virtual double mfem::DeltaCoefficient::Eval ( ElementTransformation T,
const IntegrationPoint ip 
)
inlinevirtual

A DeltaFunction cannot be evaluated. Calling this method will cause an MFEM error, terminating the application.

Implements mfem::Coefficient.

Definition at line 419 of file coefficient.hpp.

◆ EvalDelta()

double mfem::DeltaCoefficient::EvalDelta ( ElementTransformation T,
const IntegrationPoint ip 
)
virtual

The value of the function assuming we are evaluating at the delta center.

Definition at line 209 of file coefficient.cpp.

◆ GetDeltaCenter()

void mfem::DeltaCoefficient::GetDeltaCenter ( Vector center)

Write the center of the delta function into center.

Definition at line 203 of file coefficient.cpp.

◆ Scale()

double mfem::DeltaCoefficient::Scale ( )
inline

Return the scale factor times the optional time dependent function. Returns \( s T(t) \) with \( T(t) = 1 \) when not set by the user.

Definition at line 404 of file coefficient.hpp.

◆ SetDeltaCenter()

void mfem::DeltaCoefficient::SetDeltaCenter ( const Vector center)

Set the center location of the delta function.

Definition at line 195 of file coefficient.cpp.

◆ SetFunction()

void mfem::DeltaCoefficient::SetFunction ( double(*)(double)  f)
inline

Set a time-dependent function that multiplies the Scale().

Definition at line 382 of file coefficient.hpp.

◆ SetScale()

void mfem::DeltaCoefficient::SetScale ( double  s_)
inline

Set the scale value multiplying the delta function.

Definition at line 379 of file coefficient.hpp.

◆ SetTime()

void mfem::DeltaCoefficient::SetTime ( double  t)
virtual

Set the time for internally stored coefficients.

Reimplemented from mfem::Coefficient.

Definition at line 189 of file coefficient.cpp.

◆ SetTol()

void mfem::DeltaCoefficient::SetTol ( double  tol_)
inline

Set the tolerance used during projection onto GridFunction to identify the Mesh vertex where the Center() of the delta function lies. (default 1e-12)

Definition at line 387 of file coefficient.hpp.

◆ SetWeight()

void mfem::DeltaCoefficient::SetWeight ( Coefficient w)
inline

Set a weight Coefficient that multiplies the DeltaCoefficient.

The weight Coefficient multiplies the value returned by EvalDelta() but not the value returned by Scale(). The weight Coefficient is also used as the L2-weight function when projecting the DeltaCoefficient onto a GridFunction, so that the weighted integral of the projection is exactly equal to the Scale().

Definition at line 395 of file coefficient.hpp.

◆ Tol()

double mfem::DeltaCoefficient::Tol ( )
inline

Return the tolerance used to identify the mesh vertices.

Definition at line 407 of file coefficient.hpp.

◆ Weight()

Coefficient* mfem::DeltaCoefficient::Weight ( )
inline

See SetWeight() for description of the weight Coefficient.

Definition at line 410 of file coefficient.hpp.

Member Data Documentation

◆ center

double mfem::DeltaCoefficient::center[3]
protected

Definition at line 337 of file coefficient.hpp.

◆ scale

double mfem::DeltaCoefficient::scale
protected

Definition at line 337 of file coefficient.hpp.

◆ sdim

int mfem::DeltaCoefficient::sdim
protected

Definition at line 339 of file coefficient.hpp.

◆ tdf

double(* mfem::DeltaCoefficient::tdf) (double)
protected

Definition at line 340 of file coefficient.hpp.

◆ tol

double mfem::DeltaCoefficient::tol
protected

Definition at line 337 of file coefficient.hpp.

◆ weight

Coefficient* mfem::DeltaCoefficient::weight
protected

Definition at line 338 of file coefficient.hpp.


The documentation for this class was generated from the following files: