MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
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)
 
 DeltaCoefficient (real_t x, real_t s)
 Construct a delta function scaled by s and centered at (x,0.0,0.0)
 
 DeltaCoefficient (real_t x, real_t y, real_t s)
 Construct a delta function scaled by s and centered at (x,y,0.0)
 
 DeltaCoefficient (real_t x, real_t y, real_t z, real_t s)
 Construct a delta function scaled by s and centered at (x,y,z)
 
void SetTime (real_t t)
 Set the time for internally stored coefficients.
 
void SetDeltaCenter (const Vector &center)
 Set the center location of the delta function.
 
void SetScale (real_t s_)
 Set the scale value multiplying the delta function.
 
void SetFunction (real_t(*f)(real_t))
 Set a time-dependent function that multiplies the Scale().
 
void SetTol (real_t 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)
 
void SetWeight (Coefficient *w)
 Set a weight Coefficient that multiplies the DeltaCoefficient.
 
const real_tCenter ()
 
real_t 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.
 
real_t Tol ()
 Return the tolerance used to identify the mesh vertices.
 
CoefficientWeight ()
 See SetWeight() for description of the weight Coefficient.
 
void GetDeltaCenter (Vector &center)
 Write the center of the delta function into center.
 
virtual real_t EvalDelta (ElementTransformation &T, const IntegrationPoint &ip)
 The value of the function assuming we are evaluating at the delta center.
 
virtual real_t Eval (ElementTransformation &T, const IntegrationPoint &ip)
 A DeltaFunction cannot be evaluated. Calling this method will cause an MFEM error, terminating the application.
 
virtual ~DeltaCoefficient ()
 
- Public Member Functions inherited from mfem::Coefficient
 Coefficient ()
 
real_t GetTime ()
 Get the time for time dependent coefficients.
 
real_t Eval (ElementTransformation &T, const IntegrationPoint &ip, real_t t)
 Evaluate the coefficient in the element described by T at the point ip at time t.
 
virtual void Project (QuadratureFunction &qf)
 Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
 
virtual ~Coefficient ()
 

Protected Attributes

real_t center [3]
 
real_t scale
 
real_t tol
 
Coefficientweight
 
int sdim
 
real_t(* tdf )(real_t)
 
- Protected Attributes inherited from mfem::Coefficient
real_t 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 452 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 463 of file coefficient.hpp.

◆ DeltaCoefficient() [2/4]

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

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

Definition at line 470 of file coefficient.hpp.

◆ DeltaCoefficient() [3/4]

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

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

Definition at line 477 of file coefficient.hpp.

◆ DeltaCoefficient() [4/4]

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

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

Definition at line 484 of file coefficient.hpp.

◆ ~DeltaCoefficient()

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

Definition at line 539 of file coefficient.hpp.

Member Function Documentation

◆ Center()

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

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

Definition at line 517 of file coefficient.hpp.

◆ Eval()

virtual real_t 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 537 of file coefficient.hpp.

◆ EvalDelta()

real_t 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 252 of file coefficient.cpp.

◆ GetDeltaCenter()

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

Write the center of the delta function into center.

Definition at line 246 of file coefficient.cpp.

◆ Scale()

real_t 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 522 of file coefficient.hpp.

◆ SetDeltaCenter()

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

Set the center location of the delta function.

Definition at line 238 of file coefficient.cpp.

◆ SetFunction()

void mfem::DeltaCoefficient::SetFunction ( real_t(* )(real_t))
inline

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

Definition at line 500 of file coefficient.hpp.

◆ SetScale()

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

Set the scale value multiplying the delta function.

Definition at line 497 of file coefficient.hpp.

◆ SetTime()

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

Set the time for internally stored coefficients.

Reimplemented from mfem::Coefficient.

Definition at line 232 of file coefficient.cpp.

◆ SetTol()

void mfem::DeltaCoefficient::SetTol ( real_t 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 505 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 513 of file coefficient.hpp.

◆ Tol()

real_t mfem::DeltaCoefficient::Tol ( )
inline

Return the tolerance used to identify the mesh vertices.

Definition at line 525 of file coefficient.hpp.

◆ Weight()

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

See SetWeight() for description of the weight Coefficient.

Definition at line 528 of file coefficient.hpp.

Member Data Documentation

◆ center

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

Definition at line 455 of file coefficient.hpp.

◆ scale

real_t mfem::DeltaCoefficient::scale
protected

Definition at line 455 of file coefficient.hpp.

◆ sdim

int mfem::DeltaCoefficient::sdim
protected

Definition at line 457 of file coefficient.hpp.

◆ tdf

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

Definition at line 458 of file coefficient.hpp.

◆ tol

real_t mfem::DeltaCoefficient::tol
protected

Definition at line 455 of file coefficient.hpp.

◆ weight

Coefficient* mfem::DeltaCoefficient::weight
protected

Definition at line 456 of file coefficient.hpp.


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