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 "CHPFCRFFSplitVariablesAction.h" 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 : using namespace libMesh; 19 : 20 : registerMooseAction("PhaseFieldApp", CHPFCRFFSplitVariablesAction, "add_variable"); 21 : 22 : InputParameters 23 7 : CHPFCRFFSplitVariablesAction::validParams() 24 : { 25 7 : InputParameters params = Action::validParams(); 26 7 : MooseEnum familyEnum = AddVariableAction::getNonlinearVariableFamilies(); 27 7 : params.addClassDescription("Creates the L auxiliary variables, as well as a MultiApp along with " 28 : "transfers to set the variables, for the Cahn-Hilliard equation for " 29 : "the RFF form of the phase field crystal model"); 30 14 : params.addParam<MooseEnum>( 31 : "family", 32 : familyEnum, 33 : "Specifies the family of FE shape functions to use for the L variables"); 34 7 : MooseEnum orderEnum = AddVariableAction::getNonlinearVariableOrders(); 35 14 : params.addParam<MooseEnum>( 36 : "order", 37 : orderEnum, 38 : "Specifies the order of the FE shape function to use for the L variables"); 39 14 : params.addParam<Real>("scaling", 1.0, "Specifies a scaling factor to apply to the L variables"); 40 14 : params.addRequiredParam<unsigned int>( 41 : "num_L", "specifies the number of complex L variables will be solved for"); 42 14 : params.addRequiredParam<std::string>("L_name_base", "Base name for the complex L variables"); 43 14 : params.addRequiredParam<std::vector<FileName>>("sub_filenames", 44 : "This is the filename of the sub.i file"); 45 14 : params.addRequiredParam<AuxVariableName>("n_name", "Name of atomic density variable"); 46 : 47 7 : return params; 48 7 : } 49 : 50 7 : CHPFCRFFSplitVariablesAction::CHPFCRFFSplitVariablesAction(const InputParameters & params) 51 : : Action(params), 52 7 : _num_L(getParam<unsigned int>("num_L")), 53 14 : _L_name_base(getParam<std::string>("L_name_base")), 54 14 : _sub_filenames(getParam<std::vector<FileName>>("sub_filenames")), 55 14 : _n_name(getParam<AuxVariableName>("n_name")) 56 : { 57 7 : } 58 : 59 : void 60 7 : CHPFCRFFSplitVariablesAction::act() 61 : { 62 7 : ExecFlagEnum execute_options = MooseUtils::getDefaultExecFlagEnum(); 63 7 : execute_options = EXEC_TIMESTEP_BEGIN; 64 : 65 : // Setup MultiApp 66 7 : InputParameters poly_params = _factory.getValidParams("TransientMultiApp"); 67 14 : poly_params.set<MooseEnum>("app_type") = "PhaseFieldApp"; 68 7 : poly_params.set<ExecFlagEnum>("execute_on") = execute_options; 69 7 : poly_params.set<std::vector<FileName>>("input_files") = _sub_filenames; 70 7 : poly_params.set<unsigned int>("max_procs_per_app") = 1; 71 7 : poly_params.set<std::vector<Point>>("positions") = {Point()}; 72 14 : _problem->addMultiApp("TransientMultiApp", "HHEquationSolver", poly_params); 73 : 74 7 : poly_params = _factory.getValidParams("MultiAppNearestNodeTransfer"); 75 14 : poly_params.set<ExecFlagEnum>("execute_on") = execute_options; 76 21 : poly_params.set<std::vector<AuxVariableName>>("variable") = {_n_name}; 77 21 : poly_params.set<std::vector<VariableName>>("source_variable") = {_n_name}; 78 14 : poly_params.set<MultiAppName>("to_multi_app") = "HHEquationSolver"; 79 14 : _problem->addTransfer("MultiAppNearestNodeTransfer", _n_name + "_trans", poly_params); 80 : 81 : // Loop through the number of L variables 82 42 : for (unsigned int l = 0; l < _num_L; ++l) 83 : { 84 : // Create L base name 85 35 : std::string L_name = _L_name_base + Moose::stringify(l); 86 : 87 : // Create real L variable 88 35 : std::string real_name = L_name + "_real"; 89 : 90 : // Get string name for specified L variable type 91 : std::string var_type = AddVariableAction::variableType( 92 105 : FEType(Utility::string_to_enum<Order>(getParam<MooseEnum>("order")), 93 35 : Utility::string_to_enum<FEFamily>(getParam<MooseEnum>("family"))), 94 : /* is_fv = */ false, 95 35 : /* is_array = */ false); 96 : 97 : // Get params for specified L variable type 98 35 : InputParameters var_params = _factory.getValidParams(var_type); 99 : 100 35 : _problem->addAuxVariable(var_type, real_name, var_params); 101 : 102 70 : poly_params = _factory.getValidParams("MultiAppNearestNodeTransfer"); 103 105 : poly_params.set<std::vector<AuxVariableName>>("variable") = {real_name}; 104 105 : poly_params.set<std::vector<VariableName>>("source_variable") = {real_name}; 105 70 : poly_params.set<MultiAppName>("from_multi_app") = "HHEquationSolver"; 106 70 : _problem->addTransfer("MultiAppNearestNodeTransfer", real_name + "_trans", poly_params); 107 : 108 35 : if (l > 0) 109 : { 110 : // Create imaginary L variable IF l > 0 111 28 : std::string imag_name = L_name + "_imag"; 112 : 113 28 : _problem->addAuxVariable(var_type, imag_name, var_params); 114 : 115 56 : poly_params = _factory.getValidParams("MultiAppNearestNodeTransfer"); 116 84 : poly_params.set<std::vector<AuxVariableName>>("variable") = {imag_name}; 117 84 : poly_params.set<std::vector<VariableName>>("source_variable") = {imag_name}; 118 56 : poly_params.set<MultiAppName>("from_multi_app") = "HHEquationSolver"; 119 56 : _problem->addTransfer("MultiAppNearestNodeTransfer", imag_name + "_trans", poly_params); 120 : } 121 35 : } 122 14 : }