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 "Q2PAction.h"
11 :
12 : #include "Factory.h"
13 : #include "FEProblem.h"
14 : #include "Parser.h"
15 : #include "AddVariableAction.h"
16 : #include "libmesh/string_to_enum.h"
17 :
18 : registerMooseAction("RichardsApp", Q2PAction, "add_kernel");
19 :
20 : registerMooseAction("RichardsApp", Q2PAction, "add_aux_variable");
21 :
22 : registerMooseAction("RichardsApp", Q2PAction, "add_function");
23 :
24 : registerMooseAction("RichardsApp", Q2PAction, "add_postprocessor");
25 :
26 : InputParameters
27 76 : Q2PAction::validParams()
28 : {
29 152 : MooseEnum orders("CONSTANT FIRST SECOND THIRD FOURTH", "FIRST");
30 :
31 76 : InputParameters params = Action::validParams();
32 152 : params.addRequiredParam<NonlinearVariableName>("porepressure", "The porepressure variable");
33 152 : params.addRequiredParam<NonlinearVariableName>("saturation", "The water saturation variable");
34 152 : params.addRequiredParam<UserObjectName>(
35 : "water_density",
36 : "A RichardsDensity UserObject that defines the water density as a function of porepressure.");
37 152 : params.addRequiredParam<UserObjectName>(
38 : "water_relperm",
39 : "A RichardsRelPerm UserObject that defines the water relative permeability "
40 : "as a function of water saturation (eg RichardsRelPermPower).");
41 152 : params.addParam<UserObjectName>(
42 : "water_relperm_for_diffusion",
43 : "A RichardsRelPerm UserObject that defines the water relative permeability as a function of "
44 : "water saturation that will be used in the diffusivity Kernel (eg RichardsRelPermPower). If "
45 : "not given, water_relperm will be used instead, which is the most common use-case.");
46 152 : params.addRequiredParam<Real>("water_viscosity", "The water viscosity");
47 152 : params.addRequiredParam<UserObjectName>(
48 : "gas_density",
49 : "A RichardsDensity UserObject that defines the gas density as a function of porepressure.");
50 152 : params.addRequiredParam<UserObjectName>(
51 : "gas_relperm",
52 : "A RichardsRelPerm UserObject that defines the gas relative permeability as a "
53 : "function of water saturation (eg Q2PRelPermPowerGas).");
54 152 : params.addRequiredParam<Real>("gas_viscosity", "The gas viscosity");
55 152 : params.addRequiredParam<Real>("diffusivity", "The diffusivity");
56 152 : params.addParam<std::vector<OutputName>>("output_nodal_masses_to",
57 : {},
58 : "Output Nodal masses to this Output object. If you "
59 : "don't want any outputs, don't input anything here");
60 152 : params.addParam<std::vector<OutputName>>(
61 : "output_total_masses_to",
62 : {},
63 : "Output total water and gas mass to this Output object. If you "
64 : "don't want any outputs, don't input anything here");
65 152 : params.addParam<bool>("save_gas_flux_in_Q2PGasFluxResidual",
66 152 : false,
67 : "Save the residual for the "
68 : "Q2PPorepressureFlux into "
69 : "the AuxVariable called "
70 : "Q2PGasFluxResidual");
71 152 : params.addParam<bool>("save_water_flux_in_Q2PWaterFluxResidual",
72 152 : false,
73 : "Save the residual for the Q2PSaturationFlux into the AuxVariable called "
74 : "Q2PWaterFluxResidual");
75 152 : params.addParam<bool>("save_gas_Jacobian_in_Q2PGasJacobian",
76 152 : false,
77 : "Save the diagonal component of the Q2PPorepressureFlux Jacobian into the "
78 : "AuxVariable called Q2PGasJacobian");
79 152 : params.addParam<bool>("save_water_Jacobian_in_Q2PWaterJacobian",
80 152 : false,
81 : "Save the diagonal component of the Q2PSaturationFlux Jacobian into the "
82 : "AuxVariable called Q2PWaterJacobian");
83 152 : params.addParam<MooseEnum>(
84 : "ORDER",
85 : orders,
86 76 : "The order for the porepressure and saturation: " + orders.getRawNames() +
87 : " (only needed if you're calculating masses)");
88 76 : return params;
89 76 : }
90 :
91 76 : Q2PAction::Q2PAction(const InputParameters & params)
92 : : Action(params),
93 76 : _pp_var(getParam<NonlinearVariableName>("porepressure")),
94 152 : _sat_var(getParam<NonlinearVariableName>("saturation")),
95 76 : _water_density(getParam<UserObjectName>("water_density")),
96 76 : _water_relperm(getParam<UserObjectName>("water_relperm")),
97 152 : _water_relperm_for_diffusivity(isParamValid("water_relperm_for_diffusivity")
98 76 : ? getParam<UserObjectName>("water_relperm_for_diffusivity")
99 228 : : getParam<UserObjectName>("water_relperm")),
100 152 : _water_viscosity(getParam<Real>("water_viscosity")),
101 76 : _gas_density(getParam<UserObjectName>("gas_density")),
102 76 : _gas_relperm(getParam<UserObjectName>("gas_relperm")),
103 152 : _gas_viscosity(getParam<Real>("gas_viscosity")),
104 152 : _diffusivity(getParam<Real>("diffusivity")),
105 152 : _output_nodal_masses_to(getParam<std::vector<OutputName>>("output_nodal_masses_to")),
106 152 : _output_total_masses_to(getParam<std::vector<OutputName>>("output_total_masses_to")),
107 152 : _save_gas_flux_in_Q2PGasFluxResidual(getParam<bool>("save_gas_flux_in_Q2PGasFluxResidual")),
108 76 : _save_water_flux_in_Q2PWaterFluxResidual(
109 152 : getParam<bool>("save_water_flux_in_Q2PWaterFluxResidual")),
110 152 : _save_gas_Jacobian_in_Q2PGasJacobian(getParam<bool>("save_gas_Jacobian_in_Q2PGasJacobian")),
111 76 : _save_water_Jacobian_in_Q2PWaterJacobian(
112 228 : getParam<bool>("save_water_Jacobian_in_Q2PWaterJacobian"))
113 : {
114 76 : _nodal_masses_not_outputted = false;
115 76 : if (_output_nodal_masses_to.size() == 0)
116 76 : _nodal_masses_not_outputted = true;
117 :
118 76 : _total_masses_not_outputted = false;
119 76 : if (_output_total_masses_to.size() == 0)
120 24 : _total_masses_not_outputted = true;
121 :
122 76 : _no_mass_calculations = (_nodal_masses_not_outputted && _total_masses_not_outputted);
123 76 : }
124 :
125 : void
126 70 : Q2PAction::act()
127 : {
128 : // add the kernels
129 70 : if (_current_task == "add_kernel")
130 : {
131 : std::string kernel_name;
132 : std::string kernel_type;
133 28 : InputParameters params = _factory.getValidParams("Q2PNodalMass");
134 :
135 : kernel_name = "Q2P_nodal_water_mass";
136 : kernel_type = "Q2PNodalMass";
137 14 : params = _factory.getValidParams(kernel_type);
138 28 : params.set<NonlinearVariableName>("variable") = _sat_var;
139 42 : params.set<std::vector<VariableName>>("other_var") = {_pp_var};
140 14 : params.set<bool>("var_is_porepressure") = false;
141 14 : if (!_no_mass_calculations)
142 24 : params.set<std::vector<AuxVariableName>>("save_in") = {"Q2P_nodal_water_mass_divided_by_dt"};
143 14 : params.set<UserObjectName>("fluid_density") = _water_density;
144 14 : _problem->addKernel(kernel_type, kernel_name, params);
145 :
146 : kernel_name = "Q2P_nodal_gas_mass";
147 : kernel_type = "Q2PNodalMass";
148 14 : params = _factory.getValidParams(kernel_type);
149 28 : params.set<NonlinearVariableName>("variable") = _pp_var;
150 42 : params.set<std::vector<VariableName>>("other_var") = {_sat_var};
151 14 : params.set<bool>("var_is_porepressure") = true;
152 14 : if (!_no_mass_calculations)
153 24 : params.set<std::vector<AuxVariableName>>("save_in") = {"Q2P_nodal_gas_mass_divided_by_dt"};
154 14 : params.set<UserObjectName>("fluid_density") = _gas_density;
155 14 : _problem->addKernel(kernel_type, kernel_name, params);
156 :
157 : kernel_name = "Q2P_nodal_water_mass_old";
158 : kernel_type = "Q2PNegativeNodalMassOld";
159 14 : params = _factory.getValidParams(kernel_type);
160 28 : params.set<NonlinearVariableName>("variable") = _sat_var;
161 42 : params.set<std::vector<VariableName>>("other_var") = {_pp_var};
162 14 : params.set<bool>("var_is_porepressure") = false;
163 14 : params.set<UserObjectName>("fluid_density") = _water_density;
164 14 : _problem->addKernel(kernel_type, kernel_name, params);
165 :
166 : kernel_name = "Q2P_nodal_gas_mass_old";
167 : kernel_type = "Q2PNegativeNodalMassOld";
168 14 : params = _factory.getValidParams(kernel_type);
169 28 : params.set<NonlinearVariableName>("variable") = _pp_var;
170 42 : params.set<std::vector<VariableName>>("other_var") = {_sat_var};
171 14 : params.set<bool>("var_is_porepressure") = true;
172 14 : params.set<UserObjectName>("fluid_density") = _gas_density;
173 14 : _problem->addKernel(kernel_type, kernel_name, params);
174 :
175 : kernel_name = "Q2P_water_flux";
176 : kernel_type = "Q2PSaturationFlux";
177 14 : params = _factory.getValidParams(kernel_type);
178 28 : params.set<NonlinearVariableName>("variable") = _sat_var;
179 42 : params.set<std::vector<VariableName>>("porepressure_variable") = {_pp_var};
180 14 : params.set<UserObjectName>("fluid_density") = _water_density;
181 14 : params.set<UserObjectName>("fluid_relperm") = _water_relperm;
182 14 : params.set<Real>("fluid_viscosity") = _water_viscosity;
183 14 : if (_save_water_flux_in_Q2PWaterFluxResidual)
184 18 : params.set<std::vector<AuxVariableName>>("save_in") = {"Q2PWaterFluxResidual"};
185 14 : if (_save_water_Jacobian_in_Q2PWaterJacobian)
186 18 : params.set<std::vector<AuxVariableName>>("diag_save_in") = {"Q2PWaterJacobian"};
187 14 : _problem->addKernel(kernel_type, kernel_name, params);
188 :
189 : kernel_name = "Q2P_gas_flux";
190 : kernel_type = "Q2PPorepressureFlux";
191 14 : params = _factory.getValidParams(kernel_type);
192 28 : params.set<NonlinearVariableName>("variable") = _pp_var;
193 42 : params.set<std::vector<VariableName>>("saturation_variable") = {_sat_var};
194 14 : params.set<UserObjectName>("fluid_density") = _gas_density;
195 14 : params.set<UserObjectName>("fluid_relperm") = _gas_relperm;
196 14 : params.set<Real>("fluid_viscosity") = _gas_viscosity;
197 14 : if (_save_gas_flux_in_Q2PGasFluxResidual)
198 18 : params.set<std::vector<AuxVariableName>>("save_in") = {"Q2PGasFluxResidual"};
199 14 : if (_save_gas_Jacobian_in_Q2PGasJacobian)
200 18 : params.set<std::vector<AuxVariableName>>("diag_save_in") = {"Q2PGasJacobian"};
201 14 : _problem->addKernel(kernel_type, kernel_name, params);
202 :
203 : kernel_name = "Q2P_liquid_diffusion";
204 : kernel_type = "Q2PSaturationDiffusion";
205 14 : params = _factory.getValidParams(kernel_type);
206 28 : params.set<NonlinearVariableName>("variable") = _sat_var;
207 42 : params.set<std::vector<VariableName>>("porepressure_variable") = {_pp_var};
208 14 : params.set<UserObjectName>("fluid_density") = _water_density;
209 14 : params.set<UserObjectName>("fluid_relperm") = _water_relperm_for_diffusivity;
210 14 : params.set<Real>("fluid_viscosity") = _water_viscosity;
211 14 : params.set<Real>("diffusivity") = _diffusivity;
212 14 : _problem->addKernel(kernel_type, kernel_name, params);
213 14 : }
214 :
215 70 : if (_current_task == "add_aux_variable")
216 : {
217 38 : libMesh::FEType fe_type(Utility::string_to_enum<Order>(getParam<MooseEnum>("ORDER")),
218 19 : Utility::string_to_enum<libMesh::FEFamily>("LAGRANGE"));
219 19 : auto type = AddVariableAction::variableType(fe_type);
220 19 : auto var_params = _factory.getValidParams(type);
221 38 : var_params.set<MooseEnum>("family") = "LAGRANGE";
222 57 : var_params.set<MooseEnum>("order") = getParam<MooseEnum>("ORDER");
223 :
224 19 : if (!_no_mass_calculations)
225 : {
226 : // user wants nodal masses or total masses
227 13 : _problem->addAuxVariable(type, "Q2P_nodal_water_mass_divided_by_dt", var_params);
228 26 : _problem->addAuxVariable(type, "Q2P_nodal_gas_mass_divided_by_dt", var_params);
229 : }
230 19 : if (_save_gas_flux_in_Q2PGasFluxResidual)
231 22 : _problem->addAuxVariable(type, "Q2PGasFluxResidual", var_params);
232 19 : if (_save_water_flux_in_Q2PWaterFluxResidual)
233 22 : _problem->addAuxVariable(type, "Q2PWaterFluxResidual", var_params);
234 19 : if (_save_gas_Jacobian_in_Q2PGasJacobian)
235 22 : _problem->addAuxVariable(type, "Q2PGasJacobian", var_params);
236 19 : if (_save_water_Jacobian_in_Q2PWaterJacobian)
237 22 : _problem->addAuxVariable(type, "Q2PWaterJacobian", var_params);
238 19 : }
239 :
240 70 : if (_current_task == "add_function" && _output_total_masses_to.size() > 0)
241 : {
242 : // user wants total masses, so need to build Functions to do this
243 13 : InputParameters params = _factory.getValidParams("ParsedFunction");
244 :
245 26 : params.set<std::string>("value") = "a*b";
246 :
247 : std::vector<std::string> vars;
248 13 : vars.push_back("a");
249 13 : vars.push_back("b");
250 26 : params.set<std::vector<std::string>>("vars") = vars;
251 :
252 : std::vector<std::string> vals_water;
253 13 : vals_water.push_back("Q2P_mass_water_divided_by_dt");
254 13 : vals_water.push_back("Q2P_dt");
255 13 : params.set<std::vector<std::string>>("vals") = vals_water;
256 26 : _problem->addFunction("ParsedFunction", "Q2P_water_mass_fcn", params);
257 :
258 : std::vector<std::string> vals_gas;
259 13 : vals_gas.push_back("Q2P_mass_gas_divided_by_dt");
260 13 : vals_gas.push_back("Q2P_dt");
261 13 : params.set<std::vector<std::string>>("vals") = vals_gas;
262 26 : _problem->addFunction("ParsedFunction", "Q2P_gas_mass_fcn", params);
263 13 : }
264 :
265 70 : if (_current_task == "add_postprocessor" && _output_total_masses_to.size() > 0)
266 : {
267 : // user wants total masses, so need to build Postprocessors to do this
268 :
269 12 : InputParameters params = _factory.getValidParams("TimestepSize");
270 :
271 24 : params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_BEGIN;
272 36 : params.set<std::vector<OutputName>>("outputs") = {"none"};
273 24 : _problem->addPostprocessor("TimestepSize", "Q2P_dt", params);
274 :
275 24 : params = _factory.getValidParams("NodalSum");
276 36 : params.set<std::vector<OutputName>>("outputs") = {"none"};
277 36 : params.set<std::vector<VariableName>>("variable") = {"Q2P_nodal_water_mass_divided_by_dt"};
278 24 : _problem->addPostprocessor("NodalSum", "Q2P_mass_water_divided_by_dt", params);
279 :
280 24 : params = _factory.getValidParams("FunctionValuePostprocessor");
281 24 : params.set<FunctionName>("function") = "Q2P_water_mass_fcn";
282 12 : params.set<std::vector<OutputName>>("outputs") = _output_total_masses_to;
283 24 : _problem->addPostprocessor("FunctionValuePostprocessor", "mass_water", params);
284 :
285 24 : params = _factory.getValidParams("NodalSum");
286 36 : params.set<std::vector<OutputName>>("outputs") = {"none"};
287 36 : params.set<std::vector<VariableName>>("variable") = {"Q2P_nodal_gas_mass_divided_by_dt"};
288 24 : _problem->addPostprocessor("NodalSum", "Q2P_mass_gas_divided_by_dt", params);
289 :
290 24 : params = _factory.getValidParams("FunctionValuePostprocessor");
291 24 : params.set<FunctionName>("function") = "Q2P_gas_mass_fcn";
292 12 : params.set<std::vector<OutputName>>("outputs") = _output_total_masses_to;
293 24 : _problem->addPostprocessor("FunctionValuePostprocessor", "mass_gas", params);
294 12 : }
295 70 : }
|