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 "MultiSpeciesDiffusionPhysicsBase.h"
11 : #include "MatDiffusion.h"
12 : #include "MoosePreconditioner.h"
13 : #include "PetscSupport.h"
14 : #include "MooseEnumItem.h"
15 : #include "ActionComponent.h"
16 :
17 : InputParameters
18 54 : MultiSpeciesDiffusionPhysicsBase::validParams()
19 : {
20 54 : InputParameters params = PhysicsBase::validParams();
21 54 : params.addClassDescription(
22 : "Base class for creating a diffusion equation for multiple diffused species");
23 :
24 108 : params.addRequiredParam<std::vector<VariableName>>("species", "Species being diffused");
25 :
26 : // Diffusivity
27 108 : params.addParam<std::vector<MaterialPropertyName>>(
28 : "diffusivity_matprops",
29 : "Material properties defining the diffusion coefficient for each species");
30 108 : params.addParam<std::vector<MooseFunctorName>>(
31 : "diffusivity_functors", "Functors specifying the diffusivity for each species");
32 :
33 : // Source term
34 108 : params.addParam<std::vector<MooseFunctorName>>(
35 : "source_functors", "Source terms in the diffusion problem for each species");
36 108 : params.addParam<std::vector<Real>>("source_coefs", {1}, "Coefficient multiplying the source");
37 :
38 : // Boundary conditions
39 108 : params.addParam<std::vector<std::vector<BoundaryName>>>(
40 : "neumann_boundaries", {}, "Boundaries on which to apply a diffusive flux for each species");
41 108 : params.addParam<std::vector<std::vector<BoundaryName>>>(
42 : "dirichlet_boundaries", {}, "Boundaries on which to apply a fixed value for each species");
43 108 : params.addParam<std::vector<std::vector<MooseFunctorName>>>(
44 : "boundary_fluxes",
45 : {},
46 : "Functors to compute the diffusive flux on each Neumann boundary for each species");
47 108 : params.addParam<std::vector<std::vector<MooseFunctorName>>>(
48 : "boundary_values",
49 : {},
50 : "Functors to compute the diffusive flux on each Dirichlet boundary for each species");
51 108 : params.addParamNamesToGroup("neumann_boundaries dirichlet_boundaries boundary_fluxes "
52 : "boundary_values",
53 : "Boundary conditions");
54 :
55 : // Initial conditions
56 108 : params.addParam<std::vector<FunctionName>>(
57 : "initial_conditions_species", "Functions describing the initial conditions for the species");
58 :
59 : // Postprocessing
60 108 : params.addParam<std::vector<BoundaryName>>(
61 : "compute_diffusive_fluxes_on", {}, "Surfaces to compute the diffusive flux on");
62 :
63 : // Preconditioning is implemented so let's use it by default
64 108 : MooseEnum pc_options("default defer", "default");
65 108 : params.addParam<MooseEnum>(
66 : "preconditioning", pc_options, "Which preconditioning to use for this Physics");
67 :
68 108 : params.addParam<bool>(
69 : "use_automatic_differentiation",
70 108 : true,
71 : "Whether to use automatic differentiation for all the terms in the equation");
72 :
73 54 : return params;
74 54 : }
75 :
76 54 : MultiSpeciesDiffusionPhysicsBase::MultiSpeciesDiffusionPhysicsBase(
77 54 : const InputParameters & parameters)
78 : : PhysicsBase(parameters),
79 54 : _species_names(getParam<std::vector<VariableName>>("species")),
80 54 : _num_species(_species_names.size()),
81 108 : _neumann_boundaries(getParam<std::vector<std::vector<BoundaryName>>>("neumann_boundaries")),
82 108 : _dirichlet_boundaries(getParam<std::vector<std::vector<BoundaryName>>>("dirichlet_boundaries")),
83 162 : _use_ad(getParam<bool>("use_automatic_differentiation"))
84 : {
85 : // Keep track of variables
86 216 : for (const auto & var_name : _species_names)
87 : saveSolverVariableName(var_name);
88 :
89 : // Parameter checking
90 108 : checkTwoDVectorParamsSameLength<BoundaryName, MooseFunctorName>("neumann_boundaries",
91 : "boundary_fluxes");
92 108 : checkTwoDVectorParamsSameLength<BoundaryName, MooseFunctorName>("dirichlet_boundaries",
93 : "boundary_values");
94 54 : checkTwoDVectorParamsNoRespectiveOverlap<BoundaryName>(
95 : {"neumann_boundaries", "dirichlet_boundaries"});
96 108 : if (isParamSetByUser("source_coefs"))
97 108 : checkParamsBothSetOrNotSet("source_functors", "source_coefs");
98 108 : if (isParamValid("source_functors"))
99 108 : checkVectorParamsSameLength<VariableName, MooseFunctorName>("species", "source_functors");
100 108 : if (isParamValid("initial_conditions_species"))
101 36 : checkVectorParamsSameLength<VariableName, FunctionName>("species",
102 : "initial_conditions_species");
103 :
104 54 : addRequiredPhysicsTask("add_preconditioning");
105 54 : addRequiredPhysicsTask("add_ics_physics");
106 54 : }
107 :
108 : void
109 54 : MultiSpeciesDiffusionPhysicsBase::addPreconditioning()
110 : {
111 : // Use a multi-grid method, known to work for elliptic problems such as diffusion
112 54 : if (_preconditioning == "default")
113 : {
114 : // We only pass petsc options as that's all that's needed to set up the preconditioner
115 : const auto option_pair1 =
116 108 : std::make_pair<MooseEnumItem, std::string>(MooseEnumItem("-pc_type"), "hypre");
117 : const auto option_pair2 =
118 108 : std::make_pair<MooseEnumItem, std::string>(MooseEnumItem("-pc_hypre_type"), "boomeramg");
119 162 : addPetscPairsToPetscOptions({option_pair1, option_pair2});
120 54 : }
121 108 : }
122 :
123 : void
124 54 : MultiSpeciesDiffusionPhysicsBase::addPostprocessors()
125 : {
126 216 : for (const auto i : index_range(_species_names))
127 : {
128 162 : const auto & var_name = _species_names[i];
129 162 : for (const auto & boundary_name :
130 486 : getParam<std::vector<BoundaryName>>("compute_diffusive_fluxes_on"))
131 : {
132 : // Create the boundary integration of the flux
133 : const std::string pp_type =
134 243 : _use_ad ? "ADSideDiffusiveFluxIntegral" : "SideDiffusiveFluxIntegral";
135 162 : auto params = _factory.getValidParams(pp_type);
136 486 : params.set<std::vector<VariableName>>("variable") = {var_name};
137 : // FE variables require matprops for this postprocessor at the moment
138 324 : if (isParamValid("diffusivity_matprops"))
139 324 : params.set<MaterialPropertyName>("diffusivity") =
140 324 : getParam<std::vector<MaterialPropertyName>>("diffusivity_matprops")[i];
141 : // FV variables require functors
142 0 : else if (isParamValid("diffusivity_functors"))
143 0 : params.set<MooseFunctorName>("functor_diffusivity") =
144 0 : getParam<std::vector<MooseFunctorName>>("diffusivity_functors")[i];
145 : else
146 0 : mooseError("No diffusivity parameter specified");
147 486 : params.set<std::vector<BoundaryName>>("boundary") = {boundary_name};
148 : // Default to maximum computation
149 324 : params.set<ExecFlagEnum>("execute_on") = {
150 1134 : EXEC_INITIAL, EXEC_TIMESTEP_END, EXEC_NONLINEAR, EXEC_LINEAR};
151 324 : getProblem().addPostprocessor(
152 324 : pp_type, prefix() + "diffusive_flux_" + var_name + "_" + boundary_name, params);
153 162 : }
154 : }
155 216 : }
156 :
157 : void
158 0 : MultiSpeciesDiffusionPhysicsBase::addComponent(const ActionComponent & component)
159 : {
160 0 : for (const auto & block : component.blocks())
161 0 : _blocks.push_back(block);
162 0 : }
163 :
164 : void
165 54 : MultiSpeciesDiffusionPhysicsBase::addInitialConditions()
166 : {
167 54 : InputParameters params = getFactory().getValidParams("FunctionIC");
168 54 : assignBlocks(params, _blocks);
169 :
170 : // always obey the user specification of initial conditions
171 : // There are no default values, so no need to consider whether the app is restarting
172 216 : for (const auto i : index_range(_species_names))
173 : {
174 162 : const auto & var_name = _species_names[i];
175 540 : if (isParamValid("initial_conditions_species") &&
176 54 : shouldCreateIC(var_name,
177 : _blocks,
178 : /*whether IC is a default*/ false,
179 : /*error if already an IC*/ true))
180 : {
181 54 : params.set<VariableName>("variable") = var_name;
182 108 : params.set<FunctionName>("function") =
183 54 : getParam<std::vector<FunctionName>>("initial_conditions_species")[i];
184 :
185 216 : getProblem().addInitialCondition("FunctionIC", prefix() + var_name + "_ic", params);
186 : }
187 : }
188 54 : }
|