https://mooseframework.inl.gov
HHPFCRFFSplitVariablesAction.C
Go to the documentation of this file.
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 
11 #include "Factory.h"
12 #include "FEProblem.h"
13 #include "Conversion.h"
14 #include "AddVariableAction.h"
15 
16 #include "libmesh/string_to_enum.h"
17 
18 registerMooseAction("PhaseFieldApp", HHPFCRFFSplitVariablesAction, "add_variable");
19 
22 {
24  params.addClassDescription(
25  "Creates the L nonlinear variables for the Cahn-Hilliard equation for the RFF form of the "
26  "phase field crystal model, when using a split approach");
28  params.addParam<MooseEnum>(
29  "family",
30  familyEnum,
31  "Specifies the family of FE shape functions to use for the L variables");
33  params.addParam<MooseEnum>(
34  "order",
35  orderEnum,
36  "Specifies the order of the FE shape function to use for the L variables");
37  params.addParam<Real>("scaling", 1.0, "Specifies a scaling factor to apply to the L variables");
38  params.addRequiredParam<unsigned int>(
39  "num_L", "specifies the number of complex L variables will be solved for");
40  params.addRequiredParam<std::string>("L_name_base", "Base name for the complex L variables");
41  params.addParam<std::vector<SubdomainName>>(
42  "block", {}, "Block restriction for the variables and kernels");
43  return params;
44 }
45 
47  : Action(params),
48  _num_L(getParam<unsigned int>("num_L")),
49  _L_name_base(getParam<std::string>("L_name_base"))
50 {
51 }
52 
53 void
55 {
56 #ifdef DEBUG
57  Moose::err << "Inside the HHPFCRFFSplitVariablesAction Object\n";
58  Moose::err << "VariableBase: " << _L_name_base << "\torder: " << getParam<MooseEnum>("order")
59  << "\tfamily: " << getParam<MooseEnum>("family") << std::endl;
60 #endif
61 
62  // Loop through the number of L variables
63  for (unsigned int l = 0; l < _num_L; ++l)
64  {
65  // Create L base name
66  std::string L_name = _L_name_base + Moose::stringify(l);
67 
68  // Create real L variable
69  std::string real_name = L_name + "_real";
70 
71  const auto type = "MooseVariable";
72  auto params = _factory.getValidParams(type);
73  params.applySpecificParameters(_pars, {"family", "order", "block"});
74  params.set<std::vector<Real>>("scaling") = {getParam<Real>("scaling")};
75 
76  _problem->addVariable(type, real_name, params);
77 
78  if (l > 0)
79  {
80  // Create imaginary L variable IF l > 0
81  std::string imag_name = L_name + "_imag";
82 
83  _problem->addVariable(type, imag_name, params);
84  }
85  }
86 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
void applySpecificParameters(const InputParameters &common, const std::vector< std::string > &include, bool allow_private=false)
HHPFCRFFSplitVariablesAction(const InputParameters &params)
InputParameters getValidParams(const std::string &name) const
void addRequiredParam(const std::string &name, const std::string &doc_string)
static MooseEnum getNonlinearVariableFamilies()
Factory & _factory
Automatically generates all the L variables for the RFF phase field crystal model.
static InputParameters validParams()
static MooseEnum getNonlinearVariableOrders()
registerMooseAction("PhaseFieldApp", HHPFCRFFSplitVariablesAction, "add_variable")
const std::string & type() const
std::string stringify(const T &t)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const InputParameters & _pars
void addClassDescription(const std::string &doc_string)
std::shared_ptr< FEProblemBase > & _problem
void ErrorVector unsigned int