www.mooseframework.org
ConservedAction.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "ConservedAction.h"
11 // MOOSE includes
12 #include "Conversion.h"
13 #include "FEProblem.h"
14 #include "Factory.h"
15 #include "MooseObjectAction.h"
16 #include "MooseMesh.h"
17 #include "AddVariableAction.h"
18 
19 #include "libmesh/string_to_enum.h"
20 
21 registerMooseAction("PhaseFieldApp", ConservedAction, "add_variable");
22 
23 registerMooseAction("PhaseFieldApp", ConservedAction, "add_kernel");
24 
25 template <>
26 InputParameters
28 {
29  InputParameters params = validParams<Action>();
30  params.addClassDescription(
31  "Set up the variable(s) and the kernels needed for a conserved phase field variable."
32  " Note that for a direct solve, the element family and order are overwritten with hermite "
33  "and third.");
34  MooseEnum solves("DIRECT REVERSE_SPLIT FORWARD_SPLIT");
35  params.addRequiredParam<MooseEnum>("solve_type", solves, "Split or direct solve?");
36  // Get MooseEnums for the possible order/family options for this variable
37  MooseEnum families(AddVariableAction::getNonlinearVariableFamilies());
38  MooseEnum orders(AddVariableAction::getNonlinearVariableOrders());
39  params.addParam<MooseEnum>("family",
40  families,
41  "Specifies the family of FE "
42  "shape functions to use for this variable");
43  params.addParam<MooseEnum>("order",
44  orders,
45  "Specifies the order of the FE "
46  "shape function to use for this variable");
47  params.addParam<Real>("scaling", 1.0, "Specifies a scaling factor to apply to this variable");
48  params.addParam<bool>("implicit", true, "Whether kernels are implicit or not");
49  params.addParam<bool>(
50  "use_displaced_mesh", false, "Whether to use displaced mesh in the kernels");
51  params.addParamNamesToGroup("scaling implicit use_displaced_mesh", "Advanced");
52  params.addRequiredParam<MaterialPropertyName>("mobility", "The mobility used with the kernel");
53  params.addParam<std::vector<VariableName>>("args",
54  "Vector of variable arguments this kernel depends on");
55  params.addRequiredParam<MaterialPropertyName>(
56  "free_energy", "Base name of the free energy function F defined in a free energy material");
57  params.addRequiredParam<MaterialPropertyName>("kappa", "The kappa used with the kernel");
58 
59  return params;
60 }
61 
62 ConservedAction::ConservedAction(const InputParameters & params)
63  : Action(params),
64  _solve_type(getParam<MooseEnum>("solve_type").getEnum<SolveType>()),
65  _var_name(name()),
66  _scaling(getParam<Real>("scaling"))
67 {
68  switch (_solve_type)
69  {
70  case SolveType::DIRECT:
71  _fe_type = FEType(Utility::string_to_enum<Order>("THIRD"),
72  Utility::string_to_enum<FEFamily>("HERMITE"));
73  if (!parameters().isParamSetByAddParam("order") &&
74  !parameters().isParamSetByAddParam("family"))
75  mooseWarning("Order and family autoset to third and hermite in ConservedAction");
76  break;
79  _fe_type = FEType(Utility::string_to_enum<Order>(getParam<MooseEnum>("order")),
80  Utility::string_to_enum<FEFamily>(getParam<MooseEnum>("family")));
81  // Set name of chemical potential variable
82  _chempot_name = "chem_pot_" + _var_name;
83  break;
84  default:
85  paramError("solve_type", "Incorrect solve_type in ConservedAction");
86  }
87 }
88 
89 void
91 {
92  //
93  // Add variable(s)
94  //
95  if (_current_task == "add_variable")
96  {
97  auto type = AddVariableAction::determineType(_fe_type, 1);
98  auto var_params = _factory.getValidParams(type);
99  var_params.set<MooseEnum>("family") = Moose::stringify(_fe_type.family);
100  var_params.set<MooseEnum>("order") = _fe_type.order.get_order();
101  var_params.set<std::vector<Real>>("scaling") = {_scaling};
102 
103  // Create conserved variable _var_name
104  _problem->addVariable(type, _var_name, var_params);
105 
106  // Create chemical potential variable for split form
107  switch (_solve_type)
108  {
109  case SolveType::DIRECT:
110  break;
113  _problem->addVariable(type, _chempot_name, var_params);
114  }
115  }
116 
117  //
118  // Add Kernels
119  //
120  else if (_current_task == "add_kernel")
121  {
122  switch (_solve_type)
123  {
124  case SolveType::DIRECT:
125  // Add time derivative kernel
126  {
127  std::string kernel_type = "TimeDerivative";
128 
129  std::string kernel_name = _var_name + "_" + kernel_type;
130  InputParameters params = _factory.getValidParams(kernel_type);
131  params.set<NonlinearVariableName>("variable") = _var_name;
132  params.applyParameters(parameters());
133 
134  _problem->addKernel(kernel_type, kernel_name, params);
135  }
136 
137  // Add CahnHilliard kernel
138  {
139  std::string kernel_type = "CahnHilliard";
140 
141  std::string kernel_name = _var_name + "_" + kernel_type;
142  InputParameters params = _factory.getValidParams(kernel_type);
143  params.set<NonlinearVariableName>("variable") = _var_name;
144  params.set<MaterialPropertyName>("mob_name") = getParam<MaterialPropertyName>("mobility");
145  params.set<MaterialPropertyName>("f_name") =
146  getParam<MaterialPropertyName>("free_energy");
147  params.applyParameters(parameters());
148 
149  _problem->addKernel(kernel_type, kernel_name, params);
150  }
151 
152  // Add ACInterface kernel
153  {
154  std::string kernel_type = "CHInterface";
155 
156  std::string kernel_name = _var_name + "_" + kernel_type;
157  InputParameters params = _factory.getValidParams(kernel_type);
158  params.set<NonlinearVariableName>("variable") = _var_name;
159  params.set<MaterialPropertyName>("mob_name") = getParam<MaterialPropertyName>("mobility");
160  params.set<MaterialPropertyName>("kappa_name") = getParam<MaterialPropertyName>("kappa");
161  params.applyParameters(parameters());
162 
163  _problem->addKernel(kernel_type, kernel_name, params);
164  }
165  break;
166 
168  // Add time derivative kernel
169  {
170  std::string kernel_type = "CoupledTimeDerivative";
171 
172  std::string kernel_name = _var_name + "_" + kernel_type;
173  InputParameters params = _factory.getValidParams(kernel_type);
174  params.set<NonlinearVariableName>("variable") = _chempot_name;
175  params.set<std::vector<VariableName>>("v") = {_var_name};
176  params.applyParameters(parameters());
177 
178  _problem->addKernel(kernel_type, kernel_name, params);
179  }
180 
181  // Add SplitCHWRes kernel
182  {
183  std::string kernel_type = "SplitCHWRes";
184 
185  std::string kernel_name = _var_name + "_" + kernel_type;
186  InputParameters params = _factory.getValidParams(kernel_type);
187  params.set<NonlinearVariableName>("variable") = _chempot_name;
188  params.set<MaterialPropertyName>("mob_name") = getParam<MaterialPropertyName>("mobility");
189  params.applyParameters(parameters());
190 
191  _problem->addKernel(kernel_type, kernel_name, params);
192  }
193 
194  // Add SplitCHParsed kernel
195  {
196  std::string kernel_type = "SplitCHParsed";
197 
198  std::string kernel_name = _var_name + "_" + kernel_type;
199  InputParameters params = _factory.getValidParams(kernel_type);
200  params.set<NonlinearVariableName>("variable") = _var_name;
201  params.set<std::vector<VariableName>>("w") = {_chempot_name};
202  params.set<MaterialPropertyName>("f_name") =
203  getParam<MaterialPropertyName>("free_energy");
204  params.set<MaterialPropertyName>("kappa_name") = getParam<MaterialPropertyName>("kappa");
205  params.applyParameters(parameters());
206 
207  _problem->addKernel(kernel_type, kernel_name, params);
208  }
209  break;
210 
212  // Add time derivative kernel
213  {
214  std::string kernel_type = "TimeDerivative";
215 
216  std::string kernel_name = _var_name + "_" + kernel_type;
217  InputParameters params = _factory.getValidParams(kernel_type);
218  params.set<NonlinearVariableName>("variable") = _var_name;
219  params.applyParameters(parameters());
220 
221  _problem->addKernel(kernel_type, kernel_name, params);
222  }
223 
224  // Add MatDiffusion kernel for c residual
225  {
226  std::string kernel_type = "MatDiffusion";
227 
228  std::string kernel_name = _var_name + "_" + kernel_type;
229  InputParameters params = _factory.getValidParams(kernel_type);
230  params.set<NonlinearVariableName>("variable") = _var_name;
231  params.set<std::vector<VariableName>>("v") = {_chempot_name};
232  params.set<MaterialPropertyName>("diffusivity") =
233  getParam<MaterialPropertyName>("mobility");
234  params.applyParameters(parameters());
235 
236  _problem->addKernel(kernel_type, kernel_name, params);
237  }
238  // Add MatDiffusion kernel for chemical potential residual
239  {
240  std::string kernel_type = "MatDiffusion";
241 
242  std::string kernel_name = _chempot_name + "_" + kernel_type;
243  InputParameters params = _factory.getValidParams(kernel_type);
244  params.set<NonlinearVariableName>("variable") = _chempot_name;
245  params.set<std::vector<VariableName>>("v") = {_var_name};
246  params.set<MaterialPropertyName>("diffusivity") = getParam<MaterialPropertyName>("kappa");
247  params.applyParameters(parameters());
248 
249  _problem->addKernel(kernel_type, kernel_name, params);
250  }
251 
252  // Add CoupledMaterialDerivative kernel
253  {
254  std::string kernel_type = "CoupledMaterialDerivative";
255 
256  std::string kernel_name = _chempot_name + "_" + kernel_type;
257  InputParameters params = _factory.getValidParams(kernel_type);
258  params.set<NonlinearVariableName>("variable") = _chempot_name;
259  params.set<std::vector<VariableName>>("v") = {_var_name};
260  params.set<MaterialPropertyName>("f_name") =
261  getParam<MaterialPropertyName>("free_energy");
262  params.applyParameters(parameters());
263 
264  _problem->addKernel(kernel_type, kernel_name, params);
265  }
266 
267  // Add CoefReaction kernel
268  {
269  std::string kernel_type = "CoefReaction";
270 
271  std::string kernel_name = _chempot_name + "_" + kernel_type;
272  InputParameters params = _factory.getValidParams(kernel_type);
273  params.set<NonlinearVariableName>("variable") = _chempot_name;
274  params.set<Real>("coefficient") = -1.0;
275  params.applyParameters(parameters());
276 
277  _problem->addKernel(kernel_type, kernel_name, params);
278  }
279  }
280  }
281 }
ConservedAction::act
virtual void act() override
Definition: ConservedAction.C:90
ConservedAction::_var_name
const NonlinearVariableName _var_name
Name of the variable being created.
Definition: ConservedAction.h:43
registerMooseAction
registerMooseAction("PhaseFieldApp", ConservedAction, "add_variable")
ConservedAction::_solve_type
const SolveType _solve_type
Type of solve to use used in the action.
Definition: ConservedAction.h:41
ConservedAction::SolveType
SolveType
Type of solve.
Definition: ConservedAction.h:32
validParams< ConservedAction >
InputParameters validParams< ConservedAction >()
Definition: ConservedAction.C:27
ConservedAction::_chempot_name
std::string _chempot_name
Name of chemical potential variable for split solves.
Definition: ConservedAction.h:39
ConservedAction::_scaling
const Real _scaling
Scaling parameter.
Definition: ConservedAction.h:47
ConservedAction
Definition: ConservedAction.h:23
name
const std::string name
Definition: Setup.h:21
ConservedAction::_fe_type
FEType _fe_type
FEType for the variable being created.
Definition: ConservedAction.h:45
ConservedAction::SolveType::FORWARD_SPLIT
ConservedAction::SolveType::REVERSE_SPLIT
ConservedAction::SolveType::DIRECT
ConservedAction::ConservedAction
ConservedAction(const InputParameters &params)
Definition: ConservedAction.C:62
ConservedAction.h