21 #include "libmesh/system.h" 27 MultiAppGeneralFieldNearestNodeTransfer,
36 "Transfers field data at the MultiApp position by finding the value at the nearest " 37 "neighbor(s) in the origin application.");
38 params.
addParam<
unsigned int>(
"num_nearest_points",
40 "Number of nearest source (from) points will be chosen to " 41 "construct a value for the target point. All points will be " 42 "selected from the same origin mesh!");
52 params.
renameParam(
"use_nearest_app",
"assume_nearest_app_holds_nearest_location",
"");
56 std::vector<MooseEnum> source_types = {
57 MooseEnum(
"nodes centroids variable_default",
"variable_default")};
58 params.
addParam<std::vector<MooseEnum>>(
59 "source_type", source_types,
"Where to get the source values from for each source variable");
63 params.
addParam<
bool>(
"group_subapps",
65 "Whether to group source locations and values from all subapps " 66 "when considering a nearest-position algorithm");
75 _num_nearest_points(getParam<unsigned
int>(
"num_nearest_points")),
76 _group_subapps(getParam<bool>(
"group_subapps"))
80 "We do not support using both nearest positions matching and checking if target " 81 "points are within an app domain because the KDTrees for nearest-positions matching " 82 "are (currently) built with data from multiple applications.");
85 paramError(
"use_nearest_position",
"Cannot use nearest positions with mesh divisions");
91 "This option is only available for using mesh divisions or nearest positions regions");
96 "Cannot group subapps when considering nearest-location data as we would lose " 97 "track of the division index of the source locations");
101 "When using the 'nearest child application' data, the source data (positions and values) " 102 "are grouped on a per-application basis, so it cannot be agglomerated over all child " 103 "applications.\nNote that the option to use nearest applications for source restrictions, " 104 "but further split each child application's domain by regions closest to each position " 105 "(here the the child application's centroid), which could be conceived when " 106 "'group_subapps' = false, is also not available.");
115 const auto & source_types = getParam<std::vector<MooseEnum>>(
"source_type");
119 mooseError(
"Not enough source types specified for this number of variables. Source types must " 120 "be specified for transfers with multiple variables");
138 if (source_types[var_index] ==
"nodes")
140 else if (source_types[var_index] ==
"centroids")
169 mooseError(
"Source variable cannot be sampled at nodes as it is discontinuous");
187 "Transfer is projecting from nearest-nodes to centroids. This is likely causing " 188 "floating point indetermination in the results because multiple nodes are 'nearest' to " 189 "a centroid. Please consider using a ProjectionAux to build an elemental source " 190 "variable (for example constant monomial) before transferring");
194 "Transfer is projecting from nearest-centroids to nodes. This is likely causing " 195 "floating point indetermination in the results because multiple centroids are " 196 "'nearest' to a node. Please consider using a ProjectionAux to build a nodal source " 197 "variable (for example linear Lagrange) before transferring");
205 "This transfer has only been implemented with a uniform number of source mesh " 206 "divisions across all source applications");
211 const unsigned int var_index)
227 unsigned int max_leaf_size = 0;
234 for (
const auto app_i :
make_range(num_apps_per_tree))
255 for (
const auto & node : from_mesh.getMesh().local_node_ptr_range())
257 if (node->n_dofs(from_sys.
number(), from_var_num) < 1)
280 MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
283 MeshDivisionTransferUse::MATCH_SUBAPP_INDEX) &&
284 tree_division_index != node_div_index)
302 const auto dof = node->dof_number(from_sys.
number(), from_var_num, 0);
306 "Nearest-location is not implemented for this source variable type on " 307 "this mesh. Returning value at dof 0");
313 for (
auto & elem :
as_range(from_mesh.getMesh().local_elements_begin(),
314 from_mesh.getMesh().local_elements_end()))
316 if (elem->n_dofs(from_sys.
number(), from_var_num) < 1)
337 MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
339 MeshDivisionTransferUse::MATCH_SUBAPP_INDEX) &&
340 tree_division_index != elem_div_index)
344 const auto vertex_average = elem->vertex_average();
347 const auto transformed_vertex_average =
356 _local_points[i_source].push_back(transformed_vertex_average);
359 auto dof = elem->dof_number(from_sys.
number(), from_var_num, 0);
364 from_sys.
point_value(from_var_num, vertex_average, elem));
367 max_leaf_size =
std::max(max_leaf_size, from_mesh.getMaxLeafSize());
371 std::shared_ptr<KDTree> _kd_tree =
372 std::make_shared<KDTree>(
_local_points[i_source], max_leaf_size);
379 const std::vector<std::pair<Point, unsigned int>> & incoming_points,
380 std::vector<std::pair<Real, Real>> & outgoing_vals)
387 const std::vector<std::pair<Point, unsigned int>> & incoming_points,
388 std::vector<std::pair<Real, Real>> & outgoing_vals)
391 for (
const auto & [pt, mesh_div] : incoming_points)
395 bool point_found =
false;
417 Real val_sum = 0, dist_sum = 0;
418 for (
const auto index : return_index)
425 const auto new_distance = dist_sum / return_dist_sqr.size();
426 if (new_distance < outgoing_vals[i_pt].second)
427 outgoing_vals[i_pt] = {val_sum / return_index.size(), new_distance};
441 unsigned int num_equidistant_problems = 0;
459 _local_kdtrees[i_from]->neighborSearch(pt, num_search, return_index, return_dist_sqr);
460 auto num_found = return_dist_sqr.size();
469 mooseInfo(
"Search value conflict cannot find the origin point due to the " 470 "non-uniqueness of the coordinate collapsing reverse mapping");
473 std::vector<std::pair<Real, std::size_t>> zipped_nearest_points;
475 zipped_nearest_points.push_back(std::make_pair(return_dist_sqr[i], return_index[i]));
476 std::sort(zipped_nearest_points.begin(), zipped_nearest_points.end());
483 if (num_found > 1 && num_found == num_search &&
485 zipped_nearest_points[num_found - 2].first))
496 for (
const auto i :
make_range(num_search - 1))
498 auto index = zipped_nearest_points[i].second;
504 outgoing_vals[i_pt].second))
506 num_equidistant_problems++;
507 if (num_equidistant_problems > 1)
527 const Elem * elem)
const 533 const auto & node = elem->
node_ptr(i_node);
534 const auto & node_blocks =
mesh.getNodeBlockIds(*node);
535 std::set<SubdomainID> u;
536 std::set_intersection(
blocks.begin(),
540 std::inserter(u, u.begin()));
568 unsigned int kdtree_index,
unsigned int nested_loop_on_app_index)
const 575 return nested_loop_on_app_index;
615 const Point & pt)
const 630 const Point & pt,
const unsigned int mesh_div,
const unsigned int i_from)
const 638 auto position_index = i_from;
658 mooseError(
"Nearest-positions + source_app_must_contain_point not implemented");
671 "We should not be receiving point requests with an invalid " 672 "source mesh division index");
678 mesh_div != kd_div_index)
685 mesh_div != kd_div_index)
void mooseInfo(Args &&... args) const
virtual bool isNodal() const
Is this variable nodal.
registerMooseObjectRenamed("MooseApp", MultiAppGeneralFieldNearestNodeTransfer, "12/31/2024 24:00", MultiAppGeneralFieldNearestLocationTransfer)
bool inBlocks(const std::set< SubdomainID > &blocks, const MooseMesh &mesh, const Elem *elem) const override
registerMooseObject("MooseApp", MultiAppGeneralFieldNearestLocationTransfer)
void buildKDTrees(const unsigned int var_index)
const bool _group_subapps
Whether to group data when creating the nearest-point regions.
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
Function to check whether two variables are equal within an absolute tolerance.
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 ...
bool _source_app_must_contain_point
Whether the source app mesh must actually contain the points for them to be considered or whether the...
const bool _skip_coordinate_collapsing
Whether to skip coordinate collapsing (transformations of coordinates between applications using diff...
unsigned int getAppIndex(unsigned int kdtree_index, unsigned int app_index) const
Get the index of the app in the loop over the trees and the apps contributing data to each tree...
std::vector< bool > _use_zero_dof_for_value
Whether we can just use the local zero-indexed dof to get the value from the solution.
void registerConflict(unsigned int problem, dof_id_type dof_id, Point p, Real dist, bool local)
Register a potential value conflict, e.g.
unsigned int number() const
Get variable number coming from libMesh.
VariableName getFromVarName(unsigned int var_index)
Get the source variable name, with the suffix for array/vector variables.
const ExecFlagType & getCurrentExecuteOnFlag() const
Return/set the current execution flag.
virtual libMesh::FEContinuity getContinuity() const
Return the continuity of this variable.
Point getPointInLocalSourceFrame(unsigned int i_from, const Point &pt) const
Transform a point towards the local frame.
std::vector< std::unique_ptr< libMesh::PointLocatorBase > > _from_point_locators
Point locators, useful to examine point location with regards to domain restriction.
virtual libMesh::System & system()=0
Get the reference to the libMesh system.
const std::vector< VariableName > _from_var_names
Name of variables transferring from.
std::vector< bool > _source_is_nodes
Whether the source of the values is at nodes (true) or centroids (false) for each variable...
std::vector< FEProblemBase * > _to_problems
Number point_value(unsigned int var, const Point &p, const bool insist_on_success=true, const NumericVector< Number > *sol=nullptr) const
Number BetterOutOfMeshValue
FEProblemBase & _fe_problem
static InputParameters validParams()
This class provides an interface for common operations on field variables of both FE and FV types wit...
std::vector< Point > _from_positions
std::vector< const MeshDivision * > _from_mesh_divisions
Division of the origin mesh.
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
unsigned int getNumAppsPerTree() const
Number of applications which contributed nearest-locations to each KD-tree.
bool inMesh(const libMesh::PointLocatorBase *const pl, const Point &pt) const
bool closestToPosition(unsigned int pos_index, const Point &pt) const
Whether a point is closest to a position at the index specified than any other position.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
std::set< BoundaryID > _from_boundaries
Origin boundary(ies) restriction.
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...
auto max(const L &left, const R &right)
const Positions * _nearest_positions_obj
unsigned int variable_number(std::string_view var) const
unsigned int _num_sources
Number of KD-Trees to create.
An interface that allows the marking of invalid solutions during a solve.
unsigned int number() const
std::vector< MooseMesh * > _from_meshes
virtual unsigned int n_nodes() const=0
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
const std::vector< AuxVariableName > _to_var_names
Name of variables transferring to.
std::unique_ptr< NumericVector< Number > > solution
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
MultiAppGeneralFieldNearestLocationTransfer(const InputParameters ¶meters)
std::vector< std::vector< Point > > _local_points
All the nodes that meet the spatial restrictions in all the local source apps.
std::vector< std::vector< Real > > _local_values
Values of the variable being transferred at all the points in _local_points.
bool _search_value_conflicts
Whether to look for conflicts between origin points, multiple valid values for a target point...
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
std::vector< std::shared_ptr< KDTree > > _local_kdtrees
KD-Trees for all the local source apps.
unsigned int getGlobalSourceAppIndex(unsigned int i_from) const
Return the global app index from the local index in the "from-multiapp" transfer direction.
unsigned int _num_nearest_points
Number of points to consider.
unsigned int getNumPositions(bool initial=false) const
}
bool hasSecondOrderElements()
check if the mesh has SECOND order elements
static InputParameters validParams()
unsigned int INVALID_DIVISION_INDEX
Invalid subdomain id to return when outside the mesh division.
void evaluateInterpValuesNearestNode(const std::vector< std::pair< Point, unsigned int >> &incoming_points, std::vector< std::pair< Real, Real >> &outgoing_vals)
bool checkRestrictionsForSource(const Point &pt, const unsigned int valid_mesh_div, const unsigned int i_from) const
Examine all spatial restrictions that could preclude this source from being a valid source for this p...
const FEType & variable_type(const unsigned int i) const
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.
virtual void evaluateInterpValues(const std::vector< std::pair< Point, unsigned int >> &incoming_points, std::vector< std::pair< Real, Real >> &outgoing_vals) override
const Node * node_ptr(const unsigned int i) const
void mooseWarning(Args &&... args) const
Emits a warning prefixed with object name and type.
std::vector< std::unique_ptr< MultiAppCoordTransform > > _from_transforms
IntRange< T > make_range(T beg, T end)
virtual MooseMesh & mesh() override
virtual void prepareEvaluationOfInterpValues(const unsigned int var_index) override
bool _displaced_source_mesh
True if displaced mesh is used for the source mesh, otherwise false.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
std::set< SubdomainID > _from_blocks
Origin block(s) restriction.
Performs a geometric interpolation based on the values at the nearest nodes to a target location in t...
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
unsigned int getNumDivisions() const
Number of divisions (nearest-positions or source mesh divisions) used when building KD-Trees...
virtual void initialSetup() override
Method called at the beginning of the simulation for checking integrity or doing one-time setup...
SystemBase & sys()
Get the system this variable is part of.
void initialSetup() override
Method called at the beginning of the simulation for checking integrity or doing one-time setup...
bool onBoundaries(const std::set< BoundaryID > &boundaries, const MooseMesh &mesh, const Node *node) const
void computeNumSources()
Pre-compute the number of sources Number of KDTrees used to hold the locations and variable value dat...
std::vector< FEProblemBase * > _from_problems
void ErrorVector unsigned int
auto index_range(const T &sizable)
It is a general field transfer.
const MooseEnum & _from_mesh_division_behavior
How to use the origin mesh divisions to restrict the transfer.
const bool _use_nearest_app
Whether to keep track of the distance from the requested point to the app position.
const ExecFlagType EXEC_INITIAL