www.mooseframework.org
Public Member Functions | Protected Types | Protected Attributes | List of all members
ConservedAction Class Reference

#include <ConservedAction.h>

Inheritance diagram for ConservedAction:
[legend]

Public Member Functions

 ConservedAction (const InputParameters &params)
 
virtual void act () override
 

Protected Types

enum  SolveType { SolveType::DIRECT, SolveType::REVERSE_SPLIT, SolveType::FORWARD_SPLIT }
 Type of solve. More...
 

Protected Attributes

std::string _chempot_name
 Name of chemical potential variable for split solves. More...
 
const SolveType _solve_type
 Type of solve to use used in the action. More...
 
const NonlinearVariableName _var_name
 Name of the variable being created. More...
 
FEType _fe_type
 FEType for the variable being created. More...
 
const Real _scaling
 Scaling parameter. More...
 

Detailed Description

Definition at line 23 of file ConservedAction.h.

Member Enumeration Documentation

◆ SolveType

enum ConservedAction::SolveType
strongprotected

Type of solve.

Enumerator
DIRECT 
REVERSE_SPLIT 
FORWARD_SPLIT 

Definition at line 32 of file ConservedAction.h.

33  {
34  DIRECT,
35  REVERSE_SPLIT,
36  FORWARD_SPLIT
37  };

Constructor & Destructor Documentation

◆ ConservedAction()

ConservedAction::ConservedAction ( const InputParameters &  params)

Definition at line 62 of file ConservedAction.C.

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 }

Member Function Documentation

◆ act()

void ConservedAction::act ( )
overridevirtual

Definition at line 90 of file ConservedAction.C.

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 }

Member Data Documentation

◆ _chempot_name

std::string ConservedAction::_chempot_name
protected

Name of chemical potential variable for split solves.

Definition at line 39 of file ConservedAction.h.

Referenced by act(), and ConservedAction().

◆ _fe_type

FEType ConservedAction::_fe_type
protected

FEType for the variable being created.

Definition at line 45 of file ConservedAction.h.

Referenced by act(), and ConservedAction().

◆ _scaling

const Real ConservedAction::_scaling
protected

Scaling parameter.

Definition at line 47 of file ConservedAction.h.

Referenced by act().

◆ _solve_type

const SolveType ConservedAction::_solve_type
protected

Type of solve to use used in the action.

Definition at line 41 of file ConservedAction.h.

Referenced by act(), and ConservedAction().

◆ _var_name

const NonlinearVariableName ConservedAction::_var_name
protected

Name of the variable being created.

Definition at line 43 of file ConservedAction.h.

Referenced by act(), and ConservedAction().


The documentation for this class was generated from the following files:
ConservedAction::_var_name
const NonlinearVariableName _var_name
Name of the variable being created.
Definition: ConservedAction.h:43
ConservedAction::_solve_type
const SolveType _solve_type
Type of solve to use used in the action.
Definition: ConservedAction.h:41
Xfem::DIRECT
Definition: XFEM.h:41
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
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