Line data Source code
1 : //* This file is part of the MOOSE framework
2 : //* https://mooseframework.inl.gov
3 : //*
4 : //* All rights reserved, see COPYRIGHT for full restrictions
5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 : //*
7 : //* Licensed under LGPL 2.1, please see LICENSE for details
8 : //* https://www.gnu.org/licenses/lgpl-2.1.html
9 :
10 : #include "IPHDGAssemblyHelper.h"
11 : #include "MooseTypes.h"
12 : #include "MooseVariableDependencyInterface.h"
13 : #include "MooseVariableFE.h"
14 : #include "MooseVariableScalar.h"
15 : #include "SystemBase.h"
16 : #include "MooseMesh.h"
17 : #include "MooseObject.h"
18 : #include "MaterialPropertyInterface.h"
19 : #include "TransientInterface.h"
20 :
21 : using namespace libMesh;
22 :
23 : InputParameters
24 88181 : IPHDGAssemblyHelper::validParams()
25 : {
26 88181 : auto params = emptyInputParameters();
27 88181 : params.addRequiredParam<NonlinearVariableName>("face_variable", "The face variable");
28 88181 : return params;
29 0 : }
30 :
31 1303 : IPHDGAssemblyHelper::IPHDGAssemblyHelper(const MooseObject * const moose_obj,
32 : MooseVariableDependencyInterface * const mvdi,
33 : const TransientInterface * const ti,
34 : SystemBase & sys,
35 : const Assembly & assembly,
36 : const THREAD_ID tid,
37 : const std::set<SubdomainID> & block_ids,
38 1303 : const std::set<BoundaryID> & boundary_ids)
39 : : ThreeMaterialPropertyInterface(moose_obj, block_ids, boundary_ids),
40 1303 : _ti(*ti),
41 1303 : _u_var(sys.getFieldVariable<Real>(tid, moose_obj->getParam<NonlinearVariableName>("variable"))),
42 2606 : _u_face_var(sys.getFieldVariable<Real>(
43 1303 : tid, moose_obj->getParam<NonlinearVariableName>("face_variable"))),
44 1303 : _u_dof_indices(_u_var.dofIndices()),
45 1303 : _lm_u_dof_indices(_u_face_var.dofIndices()),
46 1303 : _u_sol(_u_var.adSln()),
47 1303 : _grad_u_sol(_u_var.adGradSln()),
48 1303 : _lm_u_sol(_u_face_var.adSln()),
49 1303 : _scalar_phi(_u_var.phi()),
50 1303 : _grad_scalar_phi(_u_var.gradPhi()),
51 1303 : _scalar_phi_face(_u_var.phiFace()),
52 1303 : _grad_scalar_phi_face(_u_var.gradPhiFace()),
53 1303 : _lm_phi_face(_u_face_var.phiFace()),
54 1303 : _elem_volume(assembly.elemVolume()),
55 1303 : _side_area(assembly.sideElemVolume()),
56 1303 : _ip_current_elem(assembly.elem()),
57 1303 : _ip_current_side(assembly.side()),
58 1303 : _ip_JxW(assembly.JxW()),
59 1303 : _ip_qrule(assembly.qRule()),
60 1303 : _ip_q_point(assembly.qPoints()),
61 1303 : _ip_JxW_face(assembly.JxWFace()),
62 1303 : _ip_qrule_face(assembly.qRuleFace()),
63 1303 : _ip_q_point_face(assembly.qPointsFace()),
64 2606 : _ip_normals(assembly.normals())
65 : {
66 1303 : mvdi->addMooseVariableDependency(&const_cast<MooseVariableFE<Real> &>(_u_var));
67 1303 : mvdi->addMooseVariableDependency(&const_cast<MooseVariableFE<Real> &>(_u_face_var));
68 1303 : }
69 :
70 : std::array<ADResidualsPacket, 2>
71 5067412 : IPHDGAssemblyHelper::taggingData() const
72 : {
73 5067412 : return {ADResidualsPacket{_scalar_re, _u_dof_indices, _u_var.scalingFactor()},
74 5067412 : ADResidualsPacket{_lm_re, _lm_u_dof_indices, _u_face_var.scalingFactor()}};
75 : }
76 :
77 : std::set<std::string>
78 544 : IPHDGAssemblyHelper::additionalROVariables()
79 : {
80 1088 : return {_u_face_var.name()};
81 544 : }
82 :
83 : void
84 104932 : IPHDGAssemblyHelper::lmDirichlet(const Moose::Functor<Real> & dirichlet_value)
85 : {
86 316032 : for (const auto qp : make_range(_ip_qrule_face->n_points()))
87 : {
88 633300 : const auto scalar_value = dirichlet_value(
89 211100 : Moose::ElemSideQpArg{
90 211100 : _ip_current_elem, _ip_current_side, qp, _ip_qrule_face, _ip_q_point_face[qp]},
91 211100 : _ti.determineState());
92 :
93 1774764 : for (const auto i : index_range(_lm_re))
94 1563664 : _lm_re(i) += _ip_JxW_face[qp] * (_lm_u_sol[qp] - scalar_value) * _lm_phi_face[i][qp];
95 : }
96 104932 : }
97 :
98 : void
99 10252 : IPHDGAssemblyHelper::lmPrescribedFlux(const Moose::Functor<Real> & flux_value)
100 : {
101 32016 : for (const auto qp : make_range(_ip_qrule_face->n_points()))
102 : {
103 65292 : const auto flux = flux_value(
104 21764 : Moose::ElemSideQpArg{
105 21764 : _ip_current_elem, _ip_current_side, qp, _ip_qrule_face, _ip_q_point_face[qp]},
106 21764 : _ti.determineState());
107 :
108 195876 : for (const auto i : index_range(_lm_re))
109 174112 : _lm_re(i) += _ip_JxW_face[qp] * flux * _lm_phi_face[i][qp];
110 : }
111 10252 : }
|