20 #include "libmesh/meshfree_interpolation.h" 21 #include "libmesh/numeric_vector.h" 22 #include "libmesh/system.h" 33 "Transfers the value of a variable within the master application at each sub-application " 34 "position and transfers the value to a field variable on the sub-application(s).");
36 "variable",
"The auxiliary variable to store the transferred values in.");
37 params.
addRequiredParam<VariableName>(
"source_variable",
"The variable to transfer from.");
44 _to_var_name(getParam<AuxVariableName>(
"variable")),
45 _from_var_name(getParam<VariableName>(
"source_variable"))
48 paramError(
"direction",
"This transfer is only unidirectional");
51 paramError(
"from_multi_app",
"This transfer direction has not been implemented");
63 "MultiAppVariableValueSampleTransfer");
66 "MultiAppVariableValueSampleTransfer");
73 "MultiAppVariableValueSampleTransfer::execute()", 5,
"Sampling a variable for transfer");
83 SubProblem & from_sub_problem = from_system_base.subproblem();
89 for (
unsigned int i = 0; i <
getToMultiApp()->numGlobalApps(); i++)
97 std::vector<Point> point_vec(1, multi_app_position);
100 const Elem * elem = (*pl)(multi_app_position);
107 mooseAssert(from_var.
sln().
size() == 1,
"No values in u!");
114 mooseError(
"Transfer failed to sample point value at point: ", multi_app_position);
124 unsigned int sys_num = to_sys->
number();
131 for (
const auto & node :
as_range(
mesh.localNodesBegin(),
mesh.localNodesEnd()))
133 if (node->n_dofs(sys_num, var_num) > 0)
136 dof_id_type dof = node->dof_number(sys_num, var_num, 0);
150 mooseError(
"Doesn't make sense to transfer a sampled variable's value from a MultiApp!!");
registerMooseObject("MooseApp", MultiAppVariableValueSampleTransfer)
const std::shared_ptr< MultiApp > getFromMultiApp() const
Get the MultiApp to transfer data from.
MooseEnum _current_direction
void variableIntegrityCheck(const AuxVariableName &var_name, bool is_from_multiapp) const
Utility to verify that the variable in the destination system exists.
MultiAppVariableValueSampleTransfer(const InputParameters ¶meters)
unsigned int size() const
Return the number of active items in the MultiMooseEnum.
const Parallel::Communicator & comm() const
virtual void setCurrentSubdomainID(const Elem *elem, const THREAD_ID tid)=0
const std::shared_ptr< MultiApp > getToMultiApp() const
Get the MultiApp to transfer data to.
const Parallel::Communicator & _communicator
AuxVariableName _to_var_name
Variable to sample in the parent application.
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
Base class for a system (of equations)
VariableName _from_var_name
Variable in the MultiApp to fill with the sampled value.
bool hasFromMultiApp() const
Whether the transfer owns a non-null from_multi_app.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
auto max(const L &left, const R &right)
unsigned int variable_number(std::string_view var) const
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
unsigned int number() const
unsigned int size() const
The number of elements that can currently be stored in the array.
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
void initialSetup() override
Method called at the beginning of the simulation for checking integrity or doing one-time setup...
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
static libMesh::System * find_sys(libMesh::EquationSystems &es, const std::string &var_name)
Small helper function for finding the system containing the variable.
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
virtual void initialSetup() override
Method called at the beginning of the simulation for checking integrity or doing one-time setup...
MooseVariableFieldBase & getActualFieldVariable(const THREAD_ID tid, const std::string &var_name) override
Returns the variable reference for requested MooseVariableField which may be in any system...
Samples a variable's value in the parent application domain at the point where the MultiApp (for each...
static InputParameters validParams()
void paramError(const std::string ¶m, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
static InputParameters validParams()
virtual void reinitElemPhys(const Elem *elem, const std::vector< Point > &phys_points_in_elem, const THREAD_ID tid)=0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Generic class for solving transient nonlinear problems.
virtual const FieldVariableValue & sln() const =0
MultiMooseEnum _directions
The directions this Transfer is to be executed on.
void max(const T &r, T &o, Request &req) const
Base class for all MultiAppTransfer objects.
virtual MooseMesh & mesh() override
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
virtual std::unique_ptr< libMesh::PointLocatorBase > getPointLocator() const
Proxy function to get a (sub)PointLocator from either the underlying libMesh mesh (default)...
virtual void set(const numeric_index_type i, const T value)=0
processor_id_type processor_id() const
SystemBase & sys()
Get the system this variable is part of.
processor_id_type processor_id() const
virtual void execute() override
Execute the transfer.