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 "ArrayLowerDIntegratedBC.h"
11 :
12 : #include "Assembly.h"
13 : #include "SubProblem.h"
14 : #include "SystemBase.h"
15 : #include "MooseVariableFE.h"
16 : #include "MooseVariableScalar.h"
17 :
18 : #include "libmesh/quadrature.h"
19 :
20 : InputParameters
21 28690 : ArrayLowerDIntegratedBC::validParams()
22 : {
23 28690 : InputParameters params = ArrayIntegratedBC::validParams();
24 28690 : params.addRequiredCoupledVar("lowerd_variable", "Lagrange multiplier");
25 28690 : return params;
26 0 : }
27 :
28 84 : ArrayLowerDIntegratedBC::ArrayLowerDIntegratedBC(const InputParameters & parameters)
29 : : ArrayIntegratedBC(parameters),
30 84 : _lowerd_var(*getArrayVar("lowerd_variable", 0)),
31 84 : _lambda(_is_implicit ? _lowerd_var.slnLower() : _lowerd_var.slnLowerOld()),
32 :
33 84 : _phi_lambda(_lowerd_var.phiLower()),
34 84 : _test_lambda(_lowerd_var.phiLower()),
35 168 : _work_vector(_count)
36 : {
37 84 : const auto & lower_domains = _lowerd_var.activeSubdomains();
38 168 : for (const auto & id : _mesh.boundaryLowerDBlocks())
39 84 : if (lower_domains.count(id) == 0)
40 0 : mooseDocumentedError(
41 : "moose",
42 : 29151,
43 0 : "'lowerd_variable' must be defined on the boundary lower-dimensional subdomain '" +
44 0 : _mesh.getSubdomainName(id) +
45 : "' that is added by Mesh/build_all_side_lowerd_mesh=true.\nThe check could be overly "
46 : "restrictive.");
47 :
48 168 : for (const auto & id : _var.activeSubdomains())
49 84 : if (_mesh.boundaryLowerDBlocks().count(id) > 0)
50 0 : paramError("variable",
51 0 : "Must not be defined on the boundary lower-dimensional subdomain '" +
52 0 : _mesh.getSubdomainName(id) + "'");
53 :
54 84 : if (_lowerd_var.count() != _count)
55 0 : paramError("lowerd_variable",
56 : "The number of components must be equal to the number of "
57 : "components of 'variable'");
58 84 : }
59 :
60 : void
61 5064 : ArrayLowerDIntegratedBC::computeResidual()
62 : {
63 5064 : ArrayIntegratedBC::computeResidual();
64 :
65 5064 : prepareVectorTagLower(_assembly, _lowerd_var.number());
66 :
67 25320 : for (_qp = 0; _qp < _qrule->n_points(); _qp++)
68 : {
69 20256 : initLowerDQpResidual();
70 40512 : for (_i = 0; _i < _test_lambda.size(); _i++)
71 : {
72 20256 : _work_vector.setZero();
73 20256 : computeLowerDQpResidual(_work_vector);
74 : mooseAssert(_work_vector.size() == _count,
75 : "Size of local residual is not equal to the number of array variable compoments");
76 20256 : _work_vector *= _JxW[_qp] * _coord[_qp];
77 20256 : _assembly.saveLocalArrayResidual(_local_re, _i, _test_lambda.size(), _work_vector);
78 : }
79 : }
80 :
81 5064 : accumulateTaggedLocalResidual();
82 5064 : }
83 :
84 : void
85 216 : ArrayLowerDIntegratedBC::computeJacobian()
86 : {
87 216 : ArrayIntegratedBC::computeJacobian();
88 :
89 216 : computeLowerDJacobian(Moose::LowerLower);
90 216 : }
91 :
92 : void
93 216 : ArrayLowerDIntegratedBC::computeLowerDJacobian(Moose::ConstraintJacobianType type)
94 : {
95 : mooseAssert(type == Moose::LowerLower || type == Moose::LowerPrimary ||
96 : type == Moose::PrimaryLower,
97 : "Jacobian types must have lower in computeLowerDJacobian");
98 :
99 216 : const ArrayVariableTestValue & test_space =
100 216 : (type == Moose::LowerLower || type == Moose::LowerPrimary) ? _test_lambda : _test;
101 0 : unsigned int ivar = (type == Moose::LowerLower || type == Moose::LowerPrimary)
102 216 : ? _lowerd_var.number()
103 0 : : _var.number();
104 :
105 216 : const ArrayVariableTestValue & loc_phi =
106 0 : (type == Moose::LowerLower || type == Moose::PrimaryLower) ? _phi_lambda : _phi;
107 0 : unsigned int jvar = (type == Moose::LowerLower || type == Moose::PrimaryLower)
108 216 : ? _lowerd_var.number()
109 0 : : _var.number();
110 :
111 : // need to transform the type for assembling Jacobian on boundary to be consistent with
112 : // Assembly::addJacobianLowerD() and Assembly::prepareLowerD().
113 216 : Moose::ConstraintJacobianType type_transformed =
114 : (type == Moose::LowerLower
115 216 : ? type
116 0 : : (type == Moose::LowerPrimary ? Moose::LowerSecondary : Moose::SecondaryLower));
117 216 : prepareMatrixTagLower(_assembly, ivar, jvar, type_transformed);
118 :
119 216 : if (_local_ke.n() == 0 || _local_ke.m() == 0)
120 0 : return;
121 :
122 1080 : for (_qp = 0; _qp < _qrule->n_points(); _qp++)
123 : {
124 864 : initLowerDQpJacobian(type);
125 1728 : for (_i = 0; _i < test_space.size(); _i++)
126 1728 : for (_j = 0; _j < loc_phi.size(); _j++)
127 : {
128 864 : RealEigenVector v = _JxW[_qp] * _coord[_qp] * computeLowerDQpJacobian(type);
129 864 : _assembly.saveDiagLocalArrayJacobian(
130 864 : _local_ke, _i, test_space.size(), _j, loc_phi.size(), ivar, v);
131 864 : }
132 : }
133 :
134 216 : accumulateTaggedLocalMatrix();
135 : }
136 :
137 : void
138 1188 : ArrayLowerDIntegratedBC::computeOffDiagJacobian(const unsigned int jvar_num)
139 : {
140 1188 : ArrayIntegratedBC::computeOffDiagJacobian(jvar_num);
141 :
142 1188 : computeLowerDOffDiagJacobian(Moose::LowerLower, jvar_num);
143 1188 : computeLowerDOffDiagJacobian(Moose::LowerPrimary, jvar_num);
144 1188 : computeLowerDOffDiagJacobian(Moose::PrimaryLower, jvar_num);
145 1188 : }
146 :
147 : void
148 3564 : ArrayLowerDIntegratedBC::computeLowerDOffDiagJacobian(Moose::ConstraintJacobianType type,
149 : const unsigned int jvar_num)
150 : {
151 : mooseAssert(type == Moose::LowerLower || type == Moose::LowerPrimary ||
152 : type == Moose::PrimaryLower,
153 : "Jacobian types must have lower in computeLowerDJacobian");
154 :
155 3564 : const ArrayVariableTestValue & test_space =
156 3564 : (type == Moose::LowerLower || type == Moose::LowerPrimary) ? _test_lambda : _test;
157 2376 : unsigned int ivar = (type == Moose::LowerLower || type == Moose::LowerPrimary)
158 4752 : ? _lowerd_var.number()
159 1188 : : _var.number();
160 :
161 3564 : const auto & jvar = getVariable(jvar_num);
162 :
163 : // need to transform the type for assembling Jacobian on boundary to be consistent with
164 : // Assembly::addJacobianLowerD() and Assembly::prepareLowerD().
165 3564 : Moose::ConstraintJacobianType type_transformed =
166 : (type == Moose::LowerLower
167 5940 : ? type
168 2376 : : (type == Moose::LowerPrimary ? Moose::LowerSecondary : Moose::SecondaryLower));
169 3564 : prepareMatrixTagLower(_assembly, ivar, jvar_num, type_transformed);
170 3564 : if (_local_ke.n() == 0 || _local_ke.m() == 0)
171 1728 : return;
172 :
173 1836 : if (jvar.fieldType() == Moose::VarFieldType::VAR_FIELD_STANDARD)
174 : {
175 0 : const auto & jv0 = static_cast<const MooseVariable &>(jvar);
176 : const VariableTestValue & loc_phi =
177 0 : (type == Moose::LowerLower || type == Moose::PrimaryLower) ? jv0.phiLower() : jv0.phiFace();
178 :
179 0 : for (_qp = 0; _qp < _qrule->n_points(); _qp++)
180 : {
181 0 : initLowerDQpOffDiagJacobian(type, jvar);
182 0 : for (_i = 0; _i < test_space.size(); _i++)
183 0 : for (_j = 0; _j < loc_phi.size(); _j++)
184 : {
185 0 : RealEigenMatrix v = _JxW[_qp] * _coord[_qp] * computeLowerDQpOffDiagJacobian(type, jvar);
186 0 : _assembly.saveFullLocalArrayJacobian(
187 0 : _local_ke, _i, test_space.size(), _j, loc_phi.size(), ivar, jvar_num, v);
188 0 : }
189 : }
190 : }
191 1836 : else if (jvar.fieldType() == Moose::VarFieldType::VAR_FIELD_ARRAY)
192 : {
193 1836 : const auto & jv1 = static_cast<const ArrayMooseVariable &>(jvar);
194 : const ArrayVariableTestValue & loc_phi =
195 1836 : (type == Moose::LowerLower || type == Moose::PrimaryLower) ? jv1.phiLower() : jv1.phiFace();
196 :
197 9180 : for (_qp = 0; _qp < _qrule->n_points(); _qp++)
198 : {
199 7344 : initLowerDQpOffDiagJacobian(type, jvar);
200 38016 : for (_i = 0; _i < test_space.size(); _i++)
201 80784 : for (_j = 0; _j < loc_phi.size(); _j++)
202 : {
203 50112 : RealEigenMatrix v = _JxW[_qp] * _coord[_qp] * computeLowerDQpOffDiagJacobian(type, jvar);
204 50112 : _assembly.saveFullLocalArrayJacobian(
205 50112 : _local_ke, _i, test_space.size(), _j, loc_phi.size(), ivar, jvar_num, v);
206 50112 : }
207 : }
208 : }
209 : else
210 0 : mooseError("Vector variable cannot be coupled into array DG kernel currently");
211 :
212 1836 : accumulateTaggedLocalMatrix();
213 : }
214 :
215 : RealEigenMatrix
216 40608 : ArrayLowerDIntegratedBC::computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType type,
217 : const MooseVariableFEBase & jvar)
218 : {
219 40608 : if (jvar.number() == _var.number())
220 : {
221 17280 : if (type == Moose::LowerPrimary)
222 34560 : return computeLowerDQpJacobian(type).asDiagonal();
223 : }
224 23328 : else if (jvar.number() == _lowerd_var.number())
225 : {
226 19008 : if (type == Moose::PrimaryLower || type == Moose::LowerLower)
227 38016 : return computeLowerDQpJacobian(type).asDiagonal();
228 : }
229 :
230 8640 : return RealEigenMatrix::Zero(_count, jvar.count());
231 : }
|