https://mooseframework.inl.gov
MultiAppGeneralFieldShapeEvaluationTransfer.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 
12 // MOOSE includes
13 #include "DisplacedProblem.h"
14 #include "FEProblem.h"
15 #include "MooseMesh.h"
16 #include "MooseTypes.h"
17 #include "MooseVariableFE.h"
18 #include "Positions.h"
19 
20 #include "libmesh/system.h"
21 #include "libmesh/mesh_function.h"
22 
24 
27 {
29  params.addClassDescription(
30  "Transfers field data at the MultiApp position using the finite element shape "
31  "functions from the origin application.");
32 
33  // Blanket ban on origin boundary restriction. Most shape functions have their support extend
34  // outside the boundary. For a true face variable, this parameter would make sense again
35  params.suppressParameter<std::vector<BoundaryName>>("from_boundaries");
36  // Shape function evaluations return an invalid value outside an app's domain anyway
37  params.suppressParameter<bool>("from_app_must_contain_point");
38 
39  return params;
40 }
41 
43  const InputParameters & parameters)
44  : MultiAppGeneralFieldTransfer(parameters)
45 {
46  // Nearest point isn't well defined for sending app-based data from main app to a multiapp
47  if (_nearest_positions_obj && isParamValid("to_multi_app") && !isParamValid("from_multi_app"))
48  paramError("use_nearest_position",
49  "Cannot use nearest-position algorithm when sending from the main application");
50 
51  // There's not much use for matching cell divisions when the locations have to match exactly
52  // to get a valid source variable value anyway
53  if (_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
54  _to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX)
55  paramError("from_mesh_division_usage",
56  "Matching division index is disabled for shape evaluation transfers");
57 }
58 
59 void
61  const unsigned int var_index)
62 {
63  _local_bboxes.clear();
66 
67  _local_meshfuns.clear();
69 }
70 
71 void
73  const unsigned int var_index, std::vector<libMesh::MeshFunction> & local_meshfuns)
74 {
75  local_meshfuns.reserve(_from_problems.size());
76 
77  // Construct a local mesh function for each origin problem
78  for (unsigned int i_from = 0; i_from < _from_problems.size(); ++i_from)
79  {
80  FEProblemBase & from_problem = *_from_problems[i_from];
81  MooseVariableFieldBase & from_var =
82  from_problem.getVariable(0,
83  _from_var_names[var_index],
86 
87  System & from_sys = from_var.sys().system();
88  unsigned int from_var_num = from_sys.variable_number(getFromVarName(var_index));
89 
90  local_meshfuns.emplace_back(getEquationSystem(from_problem, _displaced_source_mesh),
91  *from_sys.current_local_solution,
92  from_sys.get_dof_map(),
93  from_var_num);
94  local_meshfuns.back().init();
95  local_meshfuns.back().enable_out_of_mesh_mode(GeneralFieldTransfer::BetterOutOfMeshValue);
96  }
97 }
98 
99 void
101  const std::vector<std::pair<Point, unsigned int>> & incoming_points,
102  std::vector<std::pair<Real, Real>> & outgoing_vals)
103 {
105  _local_bboxes, _local_meshfuns, incoming_points, outgoing_vals);
106 }
107 
108 void
110  const std::vector<BoundingBox> & local_bboxes,
111  std::vector<libMesh::MeshFunction> & local_meshfuns,
112  const std::vector<std::pair<Point, unsigned int>> & incoming_points,
113  std::vector<std::pair<Real, Real>> & outgoing_vals)
114 {
115  dof_id_type i_pt = 0;
116  for (auto & [pt, mesh_div] : incoming_points)
117  {
118  bool point_found = false;
119  outgoing_vals[i_pt].second = GeneralFieldTransfer::BetterOutOfMeshValue;
120 
121  // Loop on all local origin problems until:
122  // - we've found the point in an app and the value at that point is valid
123  // - or if looking for conflicts between apps, we must check them all
124  // - or if looking for the nearest app, we also check them all
125  for (MooseIndex(_from_problems.size()) i_from = 0;
126  i_from < _from_problems.size() &&
127  (!point_found || _search_value_conflicts || _nearest_positions_obj);
128  ++i_from)
129  {
130  // Check spatial restrictions
131  Real distance = 0;
132  if (!acceptPointInOriginMesh(i_from, local_bboxes, pt, mesh_div, distance))
133  continue;
134  else
135  {
136  // Use mesh function to compute interpolation values
137  const auto from_global_num = getGlobalSourceAppIndex(i_from);
138  const auto local_pt = _from_transforms[from_global_num]->mapBack(pt);
139  auto val = (local_meshfuns[i_from])(local_pt);
140 
141  // Look for overlaps. The check is not active outside of overlap search because in that
142  // case we accept the first value from the lowest ranked process
143  // NOTE: There is no guarantee this will be the final value used among all problems
144  // but for shape evaluation we really do expect only one value to even be valid
145  if (detectConflict(val, outgoing_vals[i_pt].first, distance, outgoing_vals[i_pt].second))
146  {
147  // In the nearest-position/app mode, we save conflicts in the reference frame
149  registerConflict(i_from, 0, pt, distance, true);
150  else
151  registerConflict(i_from, 0, local_pt, distance, true);
152  }
153 
154  // No need to consider decision factors if value is invalid
156  continue;
157  else
158  point_found = true;
159 
160  // Assign value and the distance from the target point to the origin data
161  // distance only matters if multiple applications are providing a valid value
162  // which would mean they overlap in the reference space. The distance can help
163  // lift the indetermination on which value to select
164  if (distance < outgoing_vals[i_pt].second)
165  {
166  outgoing_vals[i_pt].first = val;
167  outgoing_vals[i_pt].second = distance;
168  }
169  }
170  }
171 
172  if (!point_found)
173  outgoing_vals[i_pt] = {GeneralFieldTransfer::BetterOutOfMeshValue,
175 
176  // Move to next point
177  i_pt++;
178  }
179 }
virtual void evaluateInterpValues(const std::vector< std::pair< Point, unsigned int >> &incoming_points, std::vector< std::pair< Real, Real >> &outgoing_vals) override
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:435
void registerConflict(unsigned int problem, dof_id_type dof_id, Point p, Real dist, bool local)
Register a potential value conflict, e.g.
void buildMeshFunctions(const unsigned int var_index, std::vector< libMesh::MeshFunction > &local_meshfuns)
VariableName getFromVarName(unsigned int var_index)
Get the source variable name, with the suffix for array/vector variables.
virtual libMesh::System & system()=0
Get the reference to the libMesh system.
const std::vector< VariableName > _from_var_names
Name of variables transferring from.
void evaluateInterpValuesWithMeshFunctions(const std::vector< BoundingBox > &local_bboxes, std::vector< libMesh::MeshFunction > &local_meshfuns, const std::vector< std::pair< Point, unsigned int >> &incoming_points, std::vector< std::pair< Real, Real >> &outgoing_vals)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void prepareEvaluationOfInterpValues(const unsigned int var_index) override
This class provides an interface for common operations on field variables of both FE and FV types wit...
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
Real distance(const Point &p)
virtual const MooseVariableFieldBase & getVariable(const THREAD_ID tid, const std::string &var_name, Moose::VarKindType expected_var_type=Moose::VarKindType::VAR_ANY, Moose::VarFieldType expected_var_field_type=Moose::VarFieldType::VAR_FIELD_ANY) const override
Returns the variable reference for requested variable which must be of the expected_var_type (Nonline...
unsigned int variable_number(std::string_view var) const
void suppressParameter(const std::string &name)
This method suppresses an inherited parameter so that it isn&#39;t required or valid in the derived class...
Evaluates origin shape functions to compute the target variables.
libMesh::EquationSystems & getEquationSystem(FEProblemBase &problem, bool use_displaced) const
Returns the Problem&#39;s equation system, displaced or not Be careful! If you transfer TO a displaced sy...
const bool _use_bounding_boxes
Whether to use bounding boxes to determine the applications that may receive point requests then send...
bool acceptPointInOriginMesh(unsigned int i_from, const std::vector< BoundingBox > &local_bboxes, const Point &pt, const unsigned int mesh_div, Real &distance) const
bool _search_value_conflicts
Whether to look for conflicts between origin points, multiple valid values for a target point...
unsigned int getGlobalSourceAppIndex(unsigned int i_from) const
Return the global app index from the local index in the "from-multiapp" transfer direction.
bool detectConflict(Real value_1, Real value_2, Real distance_1, Real distance_2) const
Detects whether two source values are valid and equidistant for a desired target location.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const MooseEnum & _to_mesh_division_behavior
How to use the target mesh divisions to restrict the transfer.
registerMooseObject("MooseApp", MultiAppGeneralFieldShapeEvaluationTransfer)
std::vector< std::unique_ptr< MultiAppCoordTransform > > _from_transforms
bool _displaced_source_mesh
True if displaced mesh is used for the source mesh, otherwise false.
void extractLocalFromBoundingBoxes(std::vector< BoundingBox > &local_bboxes)
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:195
SystemBase & sys()
Get the system this variable is part of.
std::vector< FEProblemBase * > _from_problems
It is a general field transfer.
uint8_t dof_id_type
const MooseEnum & _from_mesh_division_behavior
How to use the origin mesh divisions to restrict the transfer.