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 23142 : IPHDGAssemblyHelper::validParams()
25 : {
26 23142 : auto params = emptyInputParameters();
27 69426 : params.addRequiredParam<NonlinearVariableName>("face_variable", "The face variable");
28 23142 : return params;
29 0 : }
30 :
31 865 : 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 865 : const std::set<BoundaryID> & boundary_ids)
39 : : ThreeMaterialPropertyInterface(moose_obj, block_ids, boundary_ids),
40 865 : _ti(*ti),
41 865 : _u_var(sys.getFieldVariable<Real>(tid, moose_obj->getParam<NonlinearVariableName>("variable"))),
42 1730 : _u_face_var(sys.getFieldVariable<Real>(
43 1730 : tid, moose_obj->getParam<NonlinearVariableName>("face_variable"))),
44 865 : _u_dof_indices(_u_var.dofIndices()),
45 865 : _lm_u_dof_indices(_u_face_var.dofIndices()),
46 865 : _u_sol(_u_var.adSln()),
47 865 : _grad_u_sol(_u_var.adGradSln()),
48 865 : _lm_u_sol(_u_face_var.adSln()),
49 865 : _scalar_phi(_u_var.phi()),
50 865 : _grad_scalar_phi(_u_var.gradPhi()),
51 865 : _scalar_phi_face(_u_var.phiFace()),
52 865 : _grad_scalar_phi_face(_u_var.gradPhiFace()),
53 865 : _lm_phi_face(_u_face_var.phiFace()),
54 865 : _elem_volume(assembly.elemVolume()),
55 865 : _side_area(assembly.sideElemVolume()),
56 865 : _ip_current_elem(assembly.elem()),
57 865 : _ip_current_side(assembly.side()),
58 865 : _ip_current_side_elem(assembly.sideElem()),
59 865 : _ip_JxW(assembly.JxW()),
60 865 : _ip_qrule(assembly.qRule()),
61 865 : _ip_q_point(assembly.qPoints()),
62 865 : _ip_JxW_face(assembly.JxWFace()),
63 865 : _ip_qrule_face(assembly.qRuleFace()),
64 865 : _ip_q_point_face(assembly.qPointsFace()),
65 1730 : _ip_normals(assembly.normals())
66 : {
67 865 : mvdi->addMooseVariableDependency(&const_cast<MooseVariableFE<Real> &>(_u_var));
68 865 : mvdi->addMooseVariableDependency(&const_cast<MooseVariableFE<Real> &>(_u_face_var));
69 865 : }
70 :
71 : std::array<ADResidualsPacket, 2>
72 745794 : IPHDGAssemblyHelper::taggingData() const
73 : {
74 745794 : return {ADResidualsPacket{_scalar_re, _u_dof_indices, _u_var.scalingFactor()},
75 745794 : ADResidualsPacket{_lm_re, _lm_u_dof_indices, _u_face_var.scalingFactor()}};
76 : }
77 :
78 : std::set<std::string>
79 308 : IPHDGAssemblyHelper::additionalROVariables()
80 : {
81 924 : return {_u_face_var.name()};
82 308 : }
83 :
84 : void
85 29440 : IPHDGAssemblyHelper::lmDirichlet(const Moose::Functor<Real> & dirichlet_value)
86 : {
87 89168 : for (const auto qp : make_range(_ip_qrule_face->n_points()))
88 : {
89 179184 : const auto scalar_value = dirichlet_value(
90 59728 : Moose::ElemSideQpArg{
91 59728 : _ip_current_elem, _ip_current_side, qp, _ip_qrule_face, _ip_q_point_face[qp]},
92 59728 : _ti.determineState());
93 :
94 360744 : for (const auto i : index_range(_lm_re))
95 301016 : _lm_re(i) += _ip_JxW_face[qp] * (_lm_u_sol[qp] - scalar_value) * _lm_phi_face[i][qp];
96 : }
97 29440 : }
98 :
99 : void
100 9404 : IPHDGAssemblyHelper::lmPrescribedFlux(const Moose::Functor<Real> & flux_value)
101 : {
102 29112 : for (const auto qp : make_range(_ip_qrule_face->n_points()))
103 : {
104 59124 : const auto flux = flux_value(
105 19708 : Moose::ElemSideQpArg{
106 19708 : _ip_current_elem, _ip_current_side, qp, _ip_qrule_face, _ip_q_point_face[qp]},
107 19708 : _ti.determineState());
108 :
109 175804 : for (const auto i : index_range(_lm_re))
110 156096 : _lm_re(i) += _ip_JxW_face[qp] * flux * _lm_phi_face[i][qp];
111 : }
112 9404 : }
|