19 #include "libmesh/meshfree_interpolation.h" 20 #include "libmesh/system.h" 29 params.
addClassDescription(
"Transfers from a postprocessor to a scalar auxiliary variable.");
31 "from_postprocessor",
"The name of the Postprocessor to transfer the value from.");
33 "to_aux_scalar",
"The name of the scalar AuxVariable to transfer the value to.");
40 _from_pp_name(getParam<PostprocessorName>(
"from_postprocessor")),
41 _to_aux_name(getParam<VariableName>(
"to_aux_scalar"))
44 paramError(
"direction",
"This transfer is only unidirectional");
50 TIME_SECTION(
"MultiAppPostprocessorToAuxScalarTransfer::execute()",
52 "Performing transfer between a scalar variable and a postprocessor");
96 for (
unsigned int i = 0; i <
getToMultiApp()->numGlobalApps(); i++)
132 if (num_apps != scalar.
sln().size())
135 ") must be equal to the order of the scalar AuxVariable (",
164 mooseError(
"Child application allocation on parallel processes must be the same to support " 165 "siblings postprocessor to scalar variable transfer");
168 mooseError(
"Number of source and target child apps must match for siblings transfer");
virtual void execute() override
Execute the transfer.
const std::shared_ptr< MultiApp > getFromMultiApp() const
Get the MultiApp to transfer data from.
VariableName _to_aux_name
The name of the field variable to which the postprocessor is being transferred.
NumericVector< Number > & solution()
MooseEnum _current_direction
void reinit(bool reinit_for_derivative_reordering=false)
Fill out the VariableValue arrays from the system solution vector.
unsigned int size() const
Return the number of active items in the MultiMooseEnum.
const std::shared_ptr< MultiApp > getToMultiApp() const
Get the MultiApp to transfer data to.
PostprocessorName _from_pp_name
The name of the postprocessor that the transfer originates.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
MultiAppPostprocessorToAuxScalarTransfer(const InputParameters ¶meters)
Copies the value of a Postprocessor from one app to a scalar AuxVariable in another.
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
void setValues(Number value)
Set all of the values of this scalar variable to the same value.
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()
const PostprocessorValue & getPostprocessorValueByName(const PostprocessorName &name, std::size_t t_index=0) const
Get a read-only reference to the value associated with a Postprocessor that exists.
libMesh::Order order() const
Get the order of this variable Note: Order enum can be implicitly converted to unsigned int...
static InputParameters validParams()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
MultiMooseEnum _directions
The directions this Transfer is to be executed on.
Base class for all MultiAppTransfer objects.
Class for scalar variables (they are different).
IntRange< T > make_range(T beg, T end)
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
void insert(NumericVector< Number > &soln)
virtual void set(const numeric_index_type i, const Number value)=0
SystemBase & sys()
Get the system this variable is part of.
virtual void checkSiblingsTransferSupported() const override
Whether the transfer supports siblings transfer.
const VariableValue & sln() const
registerMooseObject("MooseApp", MultiAppPostprocessorToAuxScalarTransfer)