Go to the documentation of this file.
15 #include "libmesh/point.h"
16 #include "libmesh/vector_value.h"
17 #include "libmesh/elem.h"
61 virtual void compute()
override;
70 virtual T
value(
const Point & p) = 0;
99 Real
dotHelper(
const RealGradient & op1,
const RealGradient & op2) {
return op1 * op2; }
100 Real
dotHelper(
const RealTensor & op1,
const RealTensor & op2) {
return op1.contract(op2); }
104 for (
unsigned int i = 1; i < LIBMESH_DIM; ++i)
105 v += op1.col(i) * op2(i);
148 DenseVector<DataType>
_Fe;
150 DenseVector<DataType>
_Ue;
180 const std::vector<std::vector<ValueType>> *
_phi;
182 const std::vector<std::vector<GradientShapeType>> *
_dphi;
184 const std::vector<Real> *
_JxW;
std::vector< dof_id_type > _dof_indices
The global DOF indices.
Real dotHelper(const RealGradient &op1, const RealGradient &op2)
Helps perform multiplication of GradientTypes: a normal dot product for vectors and a contraction for...
InitialConditionTempl< RealEigenVector > ArrayInitialCondition
void setHermiteVertices()
set the temporary solution vector for node projections of Hermite variables
FEGenericBase< ValueType > FEBaseType
void choleskySolve(bool is_volume)
Perform the cholesky solves for edge, side, and interior projections.
DenseMatrix< Real > _Ke
Matrix storage member.
FEProblemBase & _fe_problem
dof_id_type _free_dofs
The number of free dofs.
const std::vector< std::vector< GradientShapeType > > * _dphi
pointers to shape function gradients
const InputParameters & parameters() const
Get the parameters of the object.
OutputTools< T >::OutputData DataType
InitialConditionTempl(const InputParameters ¶meters)
Constructor.
FEContinuity _cont
The type of continuity, e.g. C0, C1.
virtual void compute() override
Workhorse method for projecting the initial conditions for block initial conditions.
virtual T value(const Point &p)=0
The value of the variable at a point.
T gradientComponent(GradientType grad, unsigned int i)
OutputTools< T >::OutputShapeGradient GradientShapeType
This is a template class that implements the workhorse compute and computeNodal methods.
unsigned int _n
node counter
Real dotHelper(const RealTensor &op1, const RealTensor &op2)
std::vector< int > _free_dof
Stores the ids of the free dofs.
void choleskyAssembly(bool is_volume)
Assemble a small local system for cholesky solve.
void setCZeroVertices()
set the temporary solution vector for node projections of C0 variables
std::vector< unsigned int > _side_dofs
Side/edge DOF indices.
unsigned int _dim
the mesh dimension
MooseVariableFE< T > & _var
The variable that this initial condition is acting upon.
const std::vector< Real > * _JxW
pointers to the Jacobian * quadrature weights for current element
static InputParameters validParams()
unsigned int _qp
The current quadrature point, contains the "nth" node number when visiting nodes.
Keeps track of stuff related to assembling.
virtual GradientType gradient(const Point &)
The gradient of the variable at a point.
virtual ~InitialConditionTempl()
const std::vector< std::vector< ValueType > > * _phi
pointers to shape functions
RealEigenVector dotHelper(const RealVectorArrayValue &op1, const RealGradient &op2)
std::vector< char > _dof_is_fixed
Whether the degree of freedom is fixed (true/false)
unsigned int _n_nodes
The number of nodes for a given element.
virtual void computeNodal(const Point &p) override
Workhorse method for projecting the initial conditions for boundary restricted initial conditions.
InitialConditionTempl< Real > InitialCondition
InitialConditionTempl< RealVectorValue > VectorInitialCondition
dof_id_type _nc
number of dofs per node per variable
OutputTools< T >::OutputShape ValueType
OutputTools< T >::OutputGradient GradientType
Eigen::Matrix< Real, Eigen::Dynamic, LIBMESH_DIM > RealVectorArrayValue
const std::vector< Point > * _xyz_values
pointers to the xyz coordinates of the quadrature points for the current element
dof_id_type _current_dof
The current dof being operated on.
const Node * _current_node
The current node if the point we are evaluating at also happens to be a node.
InitialConditionBase serves as the abstract base class for InitialConditions and VectorInitialConditi...
const FEType & _fe_type
The finite element type for the IC variable.
const Moose::CoordinateSystemType & _coord_sys
The coordinate system type for this problem, references the value in Assembly.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
unsigned int _n_qp
The number of quadrature points for a given element.
DenseVector< DataType > _Ue
Linear solution vector.
DenseVector< DataType > _Fe
Linear b vector.
void setOtherCOneVertices()
set the temporary solution vector for node projections of non-Hermitian C1 variables
const Elem *const & _current_elem
The current element we are on will retrieving values at specific points in the domain.
virtual MooseVariableFEBase & variable() override
retrieves the MOOSE variable that this initial condition acts upon