16 #include "libmesh/fe_interface.h" 17 #include "libmesh/quadrature.h" 24 _fe_problem(*getCheckedPointerParam<
FEProblemBase *>(
"_fe_problem_base")),
26 _t(_fe_problem.time()),
27 _var(_sys.getActualFieldVariable<T>(parameters.
get<
THREAD_ID>(
"_tid"),
28 parameters.
get<VariableName>(
"variable"))),
31 _fe_problem.assembly(_tid, _var.kind() ==
Moose::
VAR_SOLVER ? _var.sys().number() : 0)),
32 _coord_sys(_assembly.coordSystem()),
33 _current_elem(_var.currentElem()),
34 _current_elem_volume(_assembly.elemVolume()),
35 _current_node(nullptr),
37 _fe_type(_var.feType()),
38 _dof_indices(_var.dofIndices())
59 _dim = _current_elem->dim();
61 const unsigned int n_nodes = _current_elem->n_nodes();
69 std::unique_ptr<FEBaseType> fe(FEBaseType::build(_dim, _fe_type));
72 std::unique_ptr<QBase> qrule(_fe_type.default_quadrature_rule(_dim));
73 std::unique_ptr<QBase> qedgerule(_fe_type.default_quadrature_rule(1));
74 std::unique_ptr<QBase> qsiderule(_fe_type.default_quadrature_rule(_dim - 1));
77 _phi = &fe->get_phi();
82 _cont = fe->get_continuity();
86 const std::vector<std::vector<GradientShapeType>> & ref_dphi = fe->get_dphi();
91 _JxW = &fe->get_JxW();
93 _xyz_values = &fe->get_xyz();
99 const unsigned int n_dofs = _dof_indices.size() / _var.count();
100 mooseAssert(_dof_indices.size() % _var.count() == 0,
101 "The number of degrees of freedom should be cleanly divisible by the variable count");
107 _dof_is_fixed.clear();
108 _dof_is_fixed.resize(n_dofs,
false);
110 _free_dof.resize(n_dofs, 0);
128 for (_n = 0; _n !=
n_nodes; ++_n)
130 _nc = FEInterface::n_dofs_at_node(_fe_type, _current_elem, _n);
134 auto curr_node = _current_elem->node_ptr(_n);
135 const auto & block_ids = _sys.mesh().getNodeBlockIds(*curr_node);
137 auto priority_block = *(block_ids.begin());
138 for (
auto id : block_ids)
139 if (_var.hasBlocks(
id))
145 if (!hasBlocks(priority_block) && _var.isNodal())
147 for (decltype(_nc) i = 0; i < _nc; ++i)
149 mask(_current_dof) =
false;
155 if (!_current_elem->is_vertex(_n))
165 else if (_fe_type.family ==
HERMITE)
166 setHermiteVertices();
167 else if (_cont ==
C_ONE)
168 setOtherCOneVertices();
176 _current_node =
nullptr;
178 auto & dof_map = _var.dofMap();
179 const bool add_p_level =
180 dof_map.should_p_refine(dof_map.var_group_from_var_number(_var.number()));
184 for (
unsigned int e = 0; e != _current_elem->n_edges(); ++e)
186 FEInterface::dofs_on_edge(_current_elem, _dim, _fe_type, e, _side_dofs, add_p_level);
191 for (
unsigned int i = 0; i != _side_dofs.size(); ++i)
192 if (!_dof_is_fixed[_side_dofs[i]])
193 _free_dof[_free_dofs++] = i;
200 fe->attach_quadrature_rule(qedgerule.get());
201 fe->edge_reinit(_current_elem, e);
202 _n_qp = qedgerule->n_points();
209 for (
unsigned int s = 0; s != _current_elem->n_sides(); ++s)
211 FEInterface::dofs_on_side(_current_elem, _dim, _fe_type, s, _side_dofs, add_p_level);
216 for (
unsigned int i = 0; i != _side_dofs.size(); ++i)
217 if (!_dof_is_fixed[_side_dofs[i]])
218 _free_dof[_free_dofs++] = i;
225 fe->attach_quadrature_rule(qsiderule.get());
226 fe->reinit(_current_elem, s);
227 _n_qp = qsiderule->n_points();
237 for (
unsigned int i = 0; i != n_dofs; ++i)
238 if (!_dof_is_fixed[i])
239 _free_dof[_free_dofs++] = i;
245 fe->attach_quadrature_rule(qrule.get());
246 fe->reinit(_current_elem);
247 _n_qp = qrule->n_points();
253 for (
unsigned int i = 0; i != n_dofs; ++i)
256 for (
size_t i = 0; i < mask.
size(); i++)
258 _var.setDofValue(_Ue(i), i);
261 template <
typename T>
269 _current_node = _current_elem->node_ptr(_n);
270 _Ue(_current_dof) =
value(*_current_node);
271 _dof_is_fixed[_current_dof] =
true;
280 _current_node = _current_elem->node_ptr(_n);
281 auto point_value =
value(*_current_node);
282 for (decltype(_nc) i = 0; i < _nc; ++i)
284 _Ue(_current_dof) = point_value(i);
285 _dof_is_fixed[_current_dof] =
true;
290 template <
typename T>
311 template <
typename T>
317 _current_node = _current_elem->node_ptr(_n);
318 _Ue(_current_dof) =
value(*_current_node);
319 _dof_is_fixed[_current_dof] =
true;
323 _Ue(_current_dof) = gradientComponent(grad, 0);
324 _dof_is_fixed[_current_dof] =
true;
329 Point nxminus = _current_elem->point(_n), nxplus = _current_elem->point(_n);
330 nxminus(0) -= TOLERANCE;
331 nxplus(0) += TOLERANCE;
335 _Ue(_current_dof) = gradientComponent(grad, 1);
336 _dof_is_fixed[_current_dof] =
true;
340 (gradientComponent(gxplus, 1) - gradientComponent(gxminus, 1)) / 2. / TOLERANCE;
341 _dof_is_fixed[_current_dof] =
true;
347 _Ue(_current_dof) = gradientComponent(grad, 2);
348 _dof_is_fixed[_current_dof] =
true;
352 (gradientComponent(gxplus, 2) - gradientComponent(gxminus, 2)) / 2. / TOLERANCE;
353 _dof_is_fixed[_current_dof] =
true;
356 Point nyminus = _current_elem->point(_n), nyplus = _current_elem->point(_n);
357 nyminus(1) -= TOLERANCE;
358 nyplus(1) += TOLERANCE;
363 (gradientComponent(gyplus, 2) - gradientComponent(gyminus, 2)) / 2. / TOLERANCE;
364 _dof_is_fixed[_current_dof] =
true;
367 Point nxmym = _current_elem->point(_n), nxmyp = _current_elem->point(_n),
368 nxpym = _current_elem->point(_n), nxpyp = _current_elem->point(_n);
369 nxmym(0) -= TOLERANCE;
370 nxmym(1) -= TOLERANCE;
371 nxmyp(0) -= TOLERANCE;
372 nxmyp(1) += TOLERANCE;
373 nxpym(0) += TOLERANCE;
374 nxpym(1) -= TOLERANCE;
375 nxpyp(0) += TOLERANCE;
376 nxpyp(1) += TOLERANCE;
382 (gradientComponent(gxpyp, 2) - gradientComponent(gxmyp, 2)) / 2. / TOLERANCE;
384 (gradientComponent(gxpym, 2) - gradientComponent(gxmym, 2)) / 2. / TOLERANCE;
386 _Ue(_current_dof) = (gxzplus - gxzminus) / 2. / TOLERANCE;
387 _dof_is_fixed[_current_dof] =
true;
399 template <
typename T>
407 _current_node = _current_elem->node_ptr(_n);
408 _Ue(_current_dof) =
value(*_current_node);
409 _dof_is_fixed[_current_dof] =
true;
412 for (
unsigned int i = 0; i != _dim; ++i)
414 _Ue(_current_dof) = gradientComponent(grad, i);
415 _dof_is_fixed[_current_dof] =
true;
426 template <
typename T>
431 for (_qp = 0; _qp < _n_qp; _qp++)
434 auto fineval =
value((*_xyz_values)[_qp]);
438 finegrad = gradient((*_xyz_values)[_qp]);
440 auto dofs_size = is_volume ? (_dof_indices.size() / _var.count()) : _side_dofs.size();
443 for (decltype(dofs_size) geomi = 0, freei = 0; geomi != dofs_size; ++geomi)
445 auto i = is_volume ? geomi : _side_dofs[geomi];
448 if (_dof_is_fixed[i])
450 for (decltype(dofs_size) geomj = 0, freej = 0; geomj != dofs_size; ++geomj)
452 auto j = is_volume ? geomj : _side_dofs[geomj];
453 if (_dof_is_fixed[j])
454 _Fe(freei) -= (*_phi)[i][_qp] * (*_phi)[j][_qp] * (*_JxW)[_qp] * _Ue(j);
456 _Ke(freei, freej) += (*_phi)[i][_qp] * (*_phi)[j][_qp] * (*_JxW)[_qp];
459 if (_dof_is_fixed[j])
460 _Fe(freei) -= dotHelper((*_dphi)[i][_qp], (*_dphi)[j][_qp]) * (*_JxW)[_qp] * _Ue(j);
462 _Ke(freei, freej) += dotHelper((*_dphi)[i][_qp], (*_dphi)[j][_qp]) * (*_JxW)[_qp];
464 if (!_dof_is_fixed[j])
467 _Fe(freei) += (*_phi)[i][_qp] * fineval * (*_JxW)[_qp];
469 _Fe(freei) += dotHelper(finegrad, (*_dphi)[i][_qp]) * (*_JxW)[_qp];
475 template <
typename T>
479 _Ke.resize(_free_dofs, _free_dofs);
481 _Fe.resize(_free_dofs);
484 choleskyAssembly(is_volume);
489 _Ke.cholesky_solve(_Fe, U);
492 for (
unsigned int i = 0; i != _free_dofs; ++i)
494 auto the_dof = is_volume ? _free_dof[i] : _side_dofs[_free_dof[i]];
498 _dof_is_fixed[the_dof] =
true;
506 _Ke.resize(_free_dofs, _free_dofs);
508 _Fe.resize(_free_dofs);
509 for (
unsigned int i = 0; i < _free_dofs; ++i)
510 _Fe(i).setZero(_var.count());
512 choleskyAssembly(is_volume);
517 for (
unsigned int i = 0; i < _var.count(); ++i)
520 for (
unsigned int j = 0; j < _free_dofs; ++j)
523 _Ke.cholesky_solve(v, x);
525 for (
unsigned int j = 0; j < _free_dofs; ++j)
530 for (
unsigned int i = 0; i != _free_dofs; ++i)
532 auto the_dof = is_volume ? _free_dof[i] : _side_dofs[_free_dof[i]];
534 libmesh_assert(ui.matrix().norm() < TOLERANCE || (ui - U(i)).matrix().norm() < TOLERANCE);
536 _dof_is_fixed[the_dof] =
true;
540 template <
typename T>
545 _var.computeNodalValues();
546 auto return_value =
value(p);
548 _var.setNodalValue(return_value);
552 _var.insert(_var.sys().solution());
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
OutputTools< T >::OutputData DataType
Class for stuff related to variables.
OutputTools< T >::OutputGradient GradientType
void setOtherCOneVertices()
set the temporary solution vector for node projections of non-Hermitian C1 variables ...
This is a template class that implements the workhorse compute and computeNodal methods.
void choleskyAssembly(bool is_volume)
Assemble a small local system for cholesky solve.
KOKKOS_INLINE_FUNCTION void choleskySolve(Real *const A, Real *const x, Real *const b, const unsigned int n)
Perform an in-place linear solve using Cholesky decomposition Matrix and right-hand-side vector are m...
virtual void computeNodal(const Point &p) override
Workhorse method for projecting the initial conditions for boundary restricted initial conditions...
InitialConditionBase serves as the abstract base class for InitialConditions and VectorInitialConditi...
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
const dof_id_type n_nodes
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
void setCZeroVertices()
set the temporary solution vector for node projections of C0 variables
void setHermiteVertices()
set the temporary solution vector for node projections of Hermite variables
InitialConditionTempl(const InputParameters ¶meters)
Constructor.
virtual void compute() override
Workhorse method for projecting the initial conditions for block initial conditions.
void choleskySolve(bool is_volume)
Perform the cholesky solves for edge, side, and interior projections.
virtual unsigned int size() const override final
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
T gradientComponent(GradientType grad, unsigned int i)
virtual ~InitialConditionTempl()
const Elem & get(const ElemType type_in)