libMesh
Public Member Functions | Public Attributes | Private Attributes | List of all members
libMesh::FEComputeData Class Reference

class FEComputeData hides arbitrary data to be passed to and from children of FEBase through the FEInterface::compute_data() method. More...

#include <fe_compute_data.h>

Public Member Functions

 FEComputeData (const EquationSystems &es, const Point &pin)
 Constructor. More...
 
void clear ()
 Clears the output data completely. More...
 
void init ()
 Inits the output data to default values, provided the fields are correctly resized. More...
 
void enable_derivative ()
 Enable the computation of shape gradients (dshape). More...
 
void disable_derivative ()
 Disable the computation of shape gradients (dshape). More...
 
bool need_derivative ()
 Check whether derivatives should be computed or not. More...
 

Public Attributes

const EquationSystemsequation_systems
 Const reference to the EquationSystems object that contains simulation-specific data. More...
 
const Pointp
 Holds the point where the data are to be computed. More...
 
std::vector< Numbershape
 Storage for the computed shape function values. More...
 
std::vector< Gradientdshape
 Storage for the computed shape derivative values. More...
 
std::vector< std::vector< Real > > local_transform
 Storage for local to global mapping at p. More...
 
Real phase
 Storage for the computed phase lag. More...
 
Real speed
 The wave speed. More...
 
Number frequency
 The frequency to evaluate shape functions including the wave number depending terms. More...
 

Private Attributes

bool _need_dshape
 variable indicating whether the shape-derivative should be computed or not. More...
 

Detailed Description

class FEComputeData hides arbitrary data to be passed to and from children of FEBase through the FEInterface::compute_data() method.

This enables the efficient computation of data on the finite element level, while maintaining library integrity.

Author
Daniel Dreyer
Date
2003 Helper class used with FEInterface::compute_data().

Definition at line 51 of file fe_compute_data.h.

Constructor & Destructor Documentation

◆ FEComputeData()

libMesh::FEComputeData::FEComputeData ( const EquationSystems es,
const Point pin 
)
inline

Constructor.

Takes the required input data and clears the output data using clear().

Definition at line 58 of file fe_compute_data.h.

References clear().

59  :
60  equation_systems(es),
61  p(pin),
62  _need_dshape(false)
63  {
64  this->clear();
65  }
const EquationSystems & equation_systems
Const reference to the EquationSystems object that contains simulation-specific data.
const Point & p
Holds the point where the data are to be computed.
void clear()
Clears the output data completely.
bool _need_dshape
variable indicating whether the shape-derivative should be computed or not.

Member Function Documentation

◆ clear()

void libMesh::FEComputeData::clear ( )

Clears the output data completely.

Definition at line 25 of file fe_compute_data.C.

References dshape, frequency, local_transform, phase, shape, and speed.

Referenced by FEComputeData().

26 {
27  this->shape.clear();
28  this->dshape.clear();
29  this->local_transform.clear();
30 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
31  this->phase = 0.;
32  this->speed = 1.;
33  this->frequency = 1.;
34 
35 #endif
36 }
Real phase
Storage for the computed phase lag.
std::vector< Gradient > dshape
Storage for the computed shape derivative values.
std::vector< std::vector< Real > > local_transform
Storage for local to global mapping at p.
Real speed
The wave speed.
Number frequency
The frequency to evaluate shape functions including the wave number depending terms.
std::vector< Number > shape
Storage for the computed shape function values.

◆ disable_derivative()

void libMesh::FEComputeData::disable_derivative ( )
inline

Disable the computation of shape gradients (dshape).

Default is disabled.

Definition at line 137 of file fe_compute_data.h.

References _need_dshape.

138  {_need_dshape=false; }
bool _need_dshape
variable indicating whether the shape-derivative should be computed or not.

◆ enable_derivative()

void libMesh::FEComputeData::enable_derivative ( )

Enable the computation of shape gradients (dshape).

Definition at line 71 of file fe_compute_data.C.

References _need_dshape, and dshape.

Referenced by libMesh::MeshFunction::_gradient_on_elem(), and libMesh::System::point_gradient().

72 {
73  _need_dshape=true;
74  if (!(this->dshape.empty()))
75  std::fill (this->dshape.begin(), this->dshape.end(), 0);
76 }
std::vector< Gradient > dshape
Storage for the computed shape derivative values.
bool _need_dshape
variable indicating whether the shape-derivative should be computed or not.

◆ init()

void libMesh::FEComputeData::init ( )

Inits the output data to default values, provided the fields are correctly resized.

Definition at line 40 of file fe_compute_data.C.

References equation_systems, frequency, libMesh::Parameters::get(), libMesh::Parameters::have_parameter(), libMesh::EquationSystems::parameters, phase, libMesh::Real, shape, and speed.

Referenced by libMesh::FEInterface::compute_data().

41 {
42  if (!(this->shape.empty()))
43  std::fill (this->shape.begin(), this->shape.end(), 0.);
44 
45 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
47  this->speed = this->equation_systems.parameters.get<Real>("speed");
48 
49  // ensure that the wavenumber k=2. * libMesh::pi * frequency / speed is well-defined.
50  libmesh_assert_not_equal_to(this->speed, 0);
51 
52  if (equation_systems.parameters.have_parameter<Number>("current frequency"))
53  this->frequency = this->equation_systems.parameters.get<Number>("current frequency");
54 
55 #if LIBMESH_USE_COMPLEX_NUMBERS
56  else if (equation_systems.parameters.have_parameter<Real>("current frequency"))
57  {
58  // please use the type Number instead.
59  libmesh_deprecated();
60  this->frequency = static_cast<Number> (this->equation_systems.parameters.get<Real>("current frequency"));
61  }
62 #endif
63 
64  this->phase = 0.;
65 
66 #endif //LIBMESH_ENABLE_INFINITE_ELEMENTS
67 
68 }
const EquationSystems & equation_systems
Const reference to the EquationSystems object that contains simulation-specific data.
bool have_parameter(std::string_view) const
Definition: parameters.h:395
Real phase
Storage for the computed phase lag.
Real speed
The wave speed.
const T & get(std::string_view) const
Definition: parameters.h:426
Number frequency
The frequency to evaluate shape functions including the wave number depending terms.
std::vector< Number > shape
Storage for the computed shape function values.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Parameters parameters
Data structure holding arbitrary parameters.

◆ need_derivative()

bool libMesh::FEComputeData::need_derivative ( )
inline

Check whether derivatives should be computed or not.

Definition at line 143 of file fe_compute_data.h.

References _need_dshape.

Referenced by libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), and libMesh::FEInterface::compute_data().

144  {return _need_dshape; }
bool _need_dshape
variable indicating whether the shape-derivative should be computed or not.

Member Data Documentation

◆ _need_dshape

bool libMesh::FEComputeData::_need_dshape
private

variable indicating whether the shape-derivative should be computed or not.

Default is false to save time and be compatible with elements where derivatives are not implemented/ cannot be computed.

Definition at line 152 of file fe_compute_data.h.

Referenced by disable_derivative(), enable_derivative(), and need_derivative().

◆ dshape

std::vector<Gradient> libMesh::FEComputeData::dshape

◆ equation_systems

const EquationSystems& libMesh::FEComputeData::equation_systems

Const reference to the EquationSystems object that contains simulation-specific data.

Definition at line 71 of file fe_compute_data.h.

Referenced by init().

◆ frequency

Number libMesh::FEComputeData::frequency

The frequency to evaluate shape functions including the wave number depending terms.

Use imaginary contributions for exponential damping

Definition at line 114 of file fe_compute_data.h.

Referenced by clear(), libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), and init().

◆ local_transform

std::vector<std::vector<Real> > libMesh::FEComputeData::local_transform

Storage for local to global mapping at p.

This is needed when the gradient in physical coordinates is of interest. FIXME: What kind of type should one use for it? The matrix-class don't look as if they were made for it and neither are the TensorTool-members.

Definition at line 96 of file fe_compute_data.h.

Referenced by libMesh::MeshFunction::_gradient_on_elem(), clear(), libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), libMesh::FEInterface::compute_data(), and libMesh::System::point_gradient().

◆ p

const Point& libMesh::FEComputeData::p

Holds the point where the data are to be computed.

Definition at line 76 of file fe_compute_data.h.

Referenced by libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), and libMesh::FEInterface::compute_data().

◆ phase

Real libMesh::FEComputeData::phase

Storage for the computed phase lag.

Definition at line 102 of file fe_compute_data.h.

Referenced by clear(), libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), and init().

◆ shape

std::vector<Number> libMesh::FEComputeData::shape

◆ speed

Real libMesh::FEComputeData::speed

The wave speed.

Definition at line 107 of file fe_compute_data.h.

Referenced by clear(), libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), and init().


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