https://mooseframework.inl.gov
MultiAppGeneralFieldNearestLocationTransfer.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 "FEProblem.h"
14 #include "MooseMesh.h"
15 #include "MooseTypes.h"
16 #include "MooseVariableFE.h"
17 #include "SystemBase.h"
18 #include "Positions.h"
19 #include "MooseAppCoordTransform.h"
20 
21 #include "libmesh/system.h"
22 
23 using namespace libMesh;
24 
27  MultiAppGeneralFieldNearestNodeTransfer,
28  "12/31/2024 24:00",
30 
33 {
35  params.addClassDescription(
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",
39  1,
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!");
43 
44  // Nearest node is historically more an extrapolation transfer
45  params.set<Real>("extrapolation_constant") = GeneralFieldTransfer::BetterOutOfMeshValue;
46  params.suppressParameter<Real>("extrapolation_constant");
47  // We dont keep track of both point distance to app and to nearest node, so we cannot guarantee
48  // that the nearest app (among the local apps, not globally) will hold the nearest location.
49  // However, if the user knows this is true, we can use this heuristic to reduce the number of apps
50  // that are requested to provide a candidate value. If the user is wrong, then the nearest
51  // location is used, which can be from the non-nearest app.
52  params.renameParam("use_nearest_app", "assume_nearest_app_holds_nearest_location", "");
53 
54  // the default of node/centroid switching based on the variable is causing lots of mistakes and
55  // bad results
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");
60 
61  // choose whether to include data from multiple apps when performing nearest-position/
62  // mesh-divisions based transfers
63  params.addParam<bool>("group_subapps",
64  false,
65  "Whether to group source locations and values from all subapps "
66  "when considering a nearest-position algorithm");
67 
68  return params;
69 }
70 
72  const InputParameters & parameters)
73  : MultiAppGeneralFieldTransfer(parameters),
75  _num_nearest_points(getParam<unsigned int>("num_nearest_points")),
76  _group_subapps(getParam<bool>("group_subapps"))
77 {
79  paramError("use_nearest_position",
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.");
84  (isParamValid("from_mesh_divisions") || isParamValid("to_mesh_divisions")))
85  paramError("use_nearest_position", "Cannot use nearest positions with mesh divisions");
86 
87  // Parameter checks on grouping subapp values
89  paramError(
90  "group_subapps",
91  "This option is only available for using mesh divisions or nearest positions regions");
92  else if (_group_subapps &&
93  (_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
94  _from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_SUBAPP_INDEX))
95  paramError("group_subapps",
96  "Cannot group subapps when considering nearest-location data as we would lose "
97  "track of the division index of the source locations");
98  else if (_group_subapps && _use_nearest_app)
99  paramError(
100  "group_subapps",
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.");
107 }
108 
109 void
111 {
113 
114  // Handle the source types ahead of time
115  const auto & source_types = getParam<std::vector<MooseEnum>>("source_type");
116  _source_is_nodes.resize(_from_var_names.size());
118  if (source_types.size() != _from_var_names.size())
119  mooseError("Not enough source types specified for this number of variables. Source types must "
120  "be specified for transfers with multiple variables");
121  for (const auto var_index : index_range(_from_var_names))
122  {
123  // Local app does not own any of the source problems
124  if (_from_problems.empty())
125  break;
126 
127  // Get some info on the source variable
128  FEProblemBase & from_problem = *_from_problems[0];
129  MooseVariableFieldBase & from_var =
130  from_problem.getVariable(0,
131  _from_var_names[var_index],
134  System & from_sys = from_var.sys().system();
135  const auto & fe_type = from_sys.variable_type(from_var.number());
136 
137  // Select where to get the origin values
138  if (source_types[var_index] == "nodes")
139  _source_is_nodes[var_index] = true;
140  else if (source_types[var_index] == "centroids")
141  _source_is_nodes[var_index] = false;
142  else
143  {
144  // "Nodal" variables are continuous for sure at nodes
145  if (from_var.isNodal())
146  _source_is_nodes[var_index] = true;
147  // for everyone else, we know they are continuous at centroids
148  else
149  _source_is_nodes[var_index] = false;
150  }
151 
152  // Some variables can be sampled directly at their 0 dofs
153  // - lagrange at nodes on a first order mesh
154  // - anything constant and elemental obviously has the 0-dof value at the centroid (or
155  // vertex-average). However, higher order elemental, even monomial, do not hold the centroid
156  // value at dof index 0 For example: pyramid has dof 0 at the center of the base, prism has dof
157  // 0 on an edge etc
158  if ((_source_is_nodes[var_index] && fe_type.family == LAGRANGE &&
159  !from_problem.mesh().hasSecondOrderElements()) ||
160  (!_source_is_nodes[var_index] && fe_type.order == CONSTANT))
161  _use_zero_dof_for_value[var_index] = true;
162  else
163  _use_zero_dof_for_value[var_index] = false;
164 
165  // Check with the source variable that it is not discontinuous at the source
166  if (_source_is_nodes[var_index])
167  if (from_var.getContinuity() == DISCONTINUOUS ||
168  from_var.getContinuity() == SIDE_DISCONTINUOUS)
169  mooseError("Source variable cannot be sampled at nodes as it is discontinuous");
170 
171  // Local app does not own any of the target problems
172  if (_to_problems.empty())
173  break;
174 
175  // Check with the target variable that we are not doing awful projections
176  MooseVariableFieldBase & to_var =
177  _to_problems[0]->getVariable(0,
178  _to_var_names[var_index],
181  System & to_sys = to_var.sys().system();
182  const auto & to_fe_type = to_sys.variable_type(to_var.number());
183  if (_source_is_nodes[var_index])
184  {
185  if (to_fe_type.order == CONSTANT)
186  mooseWarning(
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");
191  }
192  else if (to_var.isNodal())
193  mooseWarning(
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");
198  }
199 
200  // We need to improve the indexing if we are to allow this
201  if (!_from_mesh_divisions.empty())
202  for (const auto mesh_div : _from_mesh_divisions)
203  if (mesh_div->getNumDivisions() != _from_mesh_divisions[0]->getNumDivisions())
204  paramError("from_mesh_division",
205  "This transfer has only been implemented with a uniform number of source mesh "
206  "divisions across all source applications");
207 }
208 
209 void
211  const unsigned int var_index)
212 {
213  _local_kdtrees.clear();
214  _local_points.clear();
215  _local_values.clear();
216  buildKDTrees(var_index);
217 }
218 
219 void
221 {
223  const auto num_apps_per_tree = getNumAppsPerTree();
225  _local_points.resize(_num_sources);
226  _local_values.resize(_num_sources);
227  unsigned int max_leaf_size = 0;
228 
229  // Construct a local KDTree for each source. A source can be a single app or multiple apps
230  // combined (option for nearest-position / mesh-divisions)
231  for (const auto i_source : make_range(_num_sources))
232  {
233  // Nest a loop on apps in case multiple apps contribute to the same KD-Tree source
234  for (const auto app_i : make_range(num_apps_per_tree))
235  {
236  // Get the current app index
237  const auto i_from = getAppIndex(i_source, app_i);
238  // Current position index, if using nearest positions (not used for use_nearest_app)
239  const auto i_pos = _group_subapps ? i_source : (i_source % getNumDivisions());
240 
241  // Get access to the variable and some variable information
242  FEProblemBase & from_problem = *_from_problems[i_from];
243  auto & from_mesh = from_problem.mesh(_displaced_source_mesh);
244  MooseVariableFieldBase & from_var =
245  from_problem.getVariable(0,
246  _from_var_names[var_index],
249  System & from_sys = from_var.sys().system();
250  const unsigned int from_var_num = from_sys.variable_number(getFromVarName(var_index));
251 
252  // Build KDTree from nodes and nodal values
253  if (_source_is_nodes[var_index])
254  {
255  for (const auto & node : from_mesh.getMesh().local_node_ptr_range())
256  {
257  if (node->n_dofs(from_sys.number(), from_var_num) < 1)
258  continue;
259 
260  if (!_from_blocks.empty() && !inBlocks(_from_blocks, from_mesh, node))
261  continue;
262 
263  if (!_from_boundaries.empty() && !onBoundaries(_from_boundaries, from_mesh, node))
264  continue;
265 
266  // Handle the various source mesh divisions behaviors
267  // NOTE: This could be more efficient, as instead of rejecting points in the
268  // wrong division, we could just be adding them to the tree for the right division
269  if (!_from_mesh_divisions.empty())
270  {
271  const auto tree_division_index = i_source % getNumDivisions();
272  const auto node_div_index = _from_mesh_divisions[i_from]->divisionIndex(*node);
273 
274  // Spatial restriction is always active
275  if (node_div_index == MooseMeshDivision::INVALID_DIVISION_INDEX)
276  continue;
277  // We fill one tree per division index for matching subapp index or division index. We
278  // only accept source data from the division index
279  else if ((_from_mesh_division_behavior ==
280  MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
281  _to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
283  MeshDivisionTransferUse::MATCH_SUBAPP_INDEX) &&
284  tree_division_index != node_div_index)
285  continue;
286  }
287 
288  // Transformed node is in the reference space, as is the _nearest_positions_obj
289  const auto transformed_node = (*_from_transforms[getGlobalSourceAppIndex(i_from)])(*node);
290 
291  // Only add to the KDTree nodes that are closest to the 'position'
292  // When querying values at a target point, the KDTree associated to the closest
293  // position to the target point is queried
294  // We do not need to check the positions when using nearest app as we will assume
295  // (somewhat incorrectly) that all the points in each subapp are closer to that subapp
296  // than to any other
298  !closestToPosition(i_pos, transformed_node))
299  continue;
300 
301  _local_points[i_source].push_back(transformed_node);
302  const auto dof = node->dof_number(from_sys.number(), from_var_num, 0);
303  _local_values[i_source].push_back((*from_sys.solution)(dof));
304  if (!_use_zero_dof_for_value[var_index])
305  flagInvalidSolution(
306  "Nearest-location is not implemented for this source variable type on "
307  "this mesh. Returning value at dof 0");
308  }
309  }
310  // Build KDTree from centroids and value at centroids
311  else
312  {
313  for (auto & elem : as_range(from_mesh.getMesh().local_elements_begin(),
314  from_mesh.getMesh().local_elements_end()))
315  {
316  if (elem->n_dofs(from_sys.number(), from_var_num) < 1)
317  continue;
318 
319  if (!_from_blocks.empty() && !inBlocks(_from_blocks, elem))
320  continue;
321 
322  if (!_from_boundaries.empty() && !onBoundaries(_from_boundaries, from_mesh, elem))
323  continue;
324 
325  // Handle the various source mesh divisions behaviors
326  if (!_from_mesh_divisions.empty())
327  {
328  const auto tree_division_index = i_source % getNumDivisions();
329  const auto elem_div_index = _from_mesh_divisions[i_from]->divisionIndex(*elem);
330 
331  // Spatial restriction is always active
332  if (elem_div_index == MooseMeshDivision::INVALID_DIVISION_INDEX)
333  continue;
334  // We fill one tree per division index for matching subapp index or division index. We
335  // only accept source data from the division index
336  else if ((_from_mesh_division_behavior ==
337  MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
339  MeshDivisionTransferUse::MATCH_SUBAPP_INDEX) &&
340  tree_division_index != elem_div_index)
341  continue;
342  }
343 
344  const auto vertex_average = elem->vertex_average();
345  // Transformed element vertex average is in the reference space, as is the
346  // _nearest_positions_obj
347  const auto transformed_vertex_average =
348  (*_from_transforms[getGlobalSourceAppIndex(i_from)])(vertex_average);
349 
350  // We do not need to check the positions when using nearest app as we will assume
351  // (somewhat incorrectly) that all the points in each subapp are closer to that subapp
353  !closestToPosition(i_pos, transformed_vertex_average))
354  continue;
355 
356  _local_points[i_source].push_back(transformed_vertex_average);
357  if (_use_zero_dof_for_value[var_index])
358  {
359  auto dof = elem->dof_number(from_sys.number(), from_var_num, 0);
360  _local_values[i_source].push_back((*from_sys.solution)(dof));
361  }
362  else
363  _local_values[i_source].push_back(
364  from_sys.point_value(from_var_num, vertex_average, elem));
365  }
366  }
367  max_leaf_size = std::max(max_leaf_size, from_mesh.getMaxLeafSize());
368  }
369 
370  // Make a KDTree from the accumulated points data
371  std::shared_ptr<KDTree> _kd_tree =
372  std::make_shared<KDTree>(_local_points[i_source], max_leaf_size);
373  _local_kdtrees[i_source] = _kd_tree;
374  }
375 }
376 
377 void
379  const std::vector<std::pair<Point, unsigned int>> & incoming_points,
380  std::vector<std::pair<Real, Real>> & outgoing_vals)
381 {
382  evaluateInterpValuesNearestNode(incoming_points, outgoing_vals);
383 }
384 
385 void
387  const std::vector<std::pair<Point, unsigned int>> & incoming_points,
388  std::vector<std::pair<Real, Real>> & outgoing_vals)
389 {
390  dof_id_type i_pt = 0;
391  for (const auto & [pt, mesh_div] : incoming_points)
392  {
393  // Reset distance
394  outgoing_vals[i_pt].second = std::numeric_limits<Real>::max();
395  bool point_found = false;
396 
397  // Loop on all sources
398  for (const auto i_from : make_range(_num_sources))
399  {
400  // Examine all restrictions for the point. This source (KDTree+values) could be ruled out
401  if (!checkRestrictionsForSource(pt, mesh_div, i_from))
402  continue;
403 
404  // TODO: Pre-allocate these two work arrays. They will be regularly resized by the
405  // searches
406  std::vector<std::size_t> return_index(_num_nearest_points);
407  std::vector<Real> return_dist_sqr(_num_nearest_points);
408 
409  // KD Tree can be empty if no points are within block/boundary/bounding box restrictions
410  if (_local_kdtrees[i_from]->numberCandidatePoints())
411  {
412  point_found = true;
413  // Note that we do not need to use the transformed_pt (in the source app frame)
414  // because the KDTree has been created in the reference frame
415  _local_kdtrees[i_from]->neighborSearch(
416  pt, _num_nearest_points, return_index, return_dist_sqr);
417  Real val_sum = 0, dist_sum = 0;
418  for (const auto index : return_index)
419  {
420  val_sum += _local_values[i_from][index];
421  dist_sum += (_local_points[i_from][index] - pt).norm();
422  }
423 
424  // If the new value found is closer than for other sources, use it
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};
428  }
429  }
430 
431  // none of the source problem meshes were within the restrictions set
432  if (!point_found)
433  outgoing_vals[i_pt] = {GeneralFieldTransfer::BetterOutOfMeshValue,
435  else if (_search_value_conflicts)
436  {
437  // Local source meshes conflicts can happen two ways:
438  // - two (or more) problems' nearest nodes are equidistant to the target point
439  // - within a problem, the num_nearest_points'th closest node is as close to the target
440  // point as the num_nearest_points + 1'th closest node
441  unsigned int num_equidistant_problems = 0;
442 
443  for (const auto i_from : make_range(_num_sources))
444  {
445  // Examine all restrictions for the point. This source (KDTree+values) could be ruled out
446  if (!checkRestrictionsForSource(pt, mesh_div, i_from))
447  continue;
448 
449  // TODO: Pre-allocate these two work arrays. They will be regularly resized by the searches
450  std::vector<std::size_t> return_index(_num_nearest_points + 1);
451  std::vector<Real> return_dist_sqr(_num_nearest_points + 1);
452 
453  // NOTE: app_index is not valid if _group_subapps = true
454  const auto app_index = i_from / getNumDivisions();
455  const auto num_search = _num_nearest_points + 1;
456 
457  if (_local_kdtrees[i_from]->numberCandidatePoints())
458  {
459  _local_kdtrees[i_from]->neighborSearch(pt, num_search, return_index, return_dist_sqr);
460  auto num_found = return_dist_sqr.size();
461 
462  // Local coordinates only accessible when not using nearest-position
463  // as we did not keep the index of the source app, only the position index
464  const Point local_pt = getPointInLocalSourceFrame(app_index, pt);
465 
466  if (!_nearest_positions_obj &&
467  _from_transforms[getGlobalSourceAppIndex(app_index)]->hasCoordinateSystemTypeChange())
469  mooseInfo("Search value conflict cannot find the origin point due to the "
470  "non-uniqueness of the coordinate collapsing reverse mapping");
471 
472  // Look for too many equidistant nodes within a problem. First zip then sort by distance
473  std::vector<std::pair<Real, std::size_t>> zipped_nearest_points;
474  for (const auto i : make_range(num_found))
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());
477 
478  // If two furthest are equally far from target point, then we have an indetermination in
479  // what is sent in this communication round from this process. However, it may not
480  // materialize to an actual conflict, as values sent from another process for the
481  // desired target point could be closer (nearest). There is no way to know at this point
482  // in the communication that a closer value exists somewhere else
483  if (num_found > 1 && num_found == num_search &&
484  MooseUtils::absoluteFuzzyEqual(zipped_nearest_points[num_found - 1].first,
485  zipped_nearest_points[num_found - 2].first))
486  {
488  registerConflict(app_index, 0, pt, outgoing_vals[i_pt].second, true);
489  else
490  registerConflict(app_index, 0, local_pt, outgoing_vals[i_pt].second, true);
491  }
492 
493  // Recompute the distance for this problem. If it matches the cached value more than
494  // once it means multiple problems provide equidistant values for this point
495  Real dist_sum = 0;
496  for (const auto i : make_range(num_search - 1))
497  {
498  auto index = zipped_nearest_points[i].second;
499  dist_sum += (_local_points[i_from][index] - pt).norm();
500  }
501 
502  // Compare to the selected value found after looking at all the problems
503  if (MooseUtils::absoluteFuzzyEqual(dist_sum / return_dist_sqr.size(),
504  outgoing_vals[i_pt].second))
505  {
506  num_equidistant_problems++;
507  if (num_equidistant_problems > 1)
508  {
510  registerConflict(app_index, 0, pt, outgoing_vals[i_pt].second, true);
511  else
512  registerConflict(app_index, 0, local_pt, outgoing_vals[i_pt].second, true);
513  }
514  }
515  }
516  }
517  }
518 
519  // Move to next point
520  i_pt++;
521  }
522 }
523 
524 bool
526  const MooseMesh & mesh,
527  const Elem * elem) const
528 {
529  // We need to override the definition of block restriction for an element
530  // because we have to consider whether each node of an element is adjacent to a block
531  for (const auto & i_node : make_range(elem->n_nodes()))
532  {
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(),
537  blocks.end(),
538  node_blocks.begin(),
539  node_blocks.end(),
540  std::inserter(u, u.begin()));
541  if (!u.empty())
542  return true;
543  }
544  return false;
545 }
546 
547 void
549 {
550  // Number of source = number of KDTrees.
551  // Using mesh divisions or nearest-positions, for every app we use 1 tree per division
552  if (!_from_mesh_divisions.empty() ||
555  // If we group apps, then we only use one tree per division (nearest-position region)
558  EXEC_INITIAL);
559  // Regular case: 1 KDTree per app
560  // Also if use_nearest_app = true, the number of problems is better than the number of positions,
561  // because some of the positions are positions of child applications that are not local
562  else
563  _num_sources = _from_problems.size();
564 }
565 
566 unsigned int
568  unsigned int kdtree_index, unsigned int nested_loop_on_app_index) const
569 {
570  // Each app is mapped to a single KD Tree
571  if (_use_nearest_app)
572  return kdtree_index;
573  // We are looping over all the apps that are grouped together
574  else if (_group_subapps)
575  return nested_loop_on_app_index;
576  // There are num_divisions trees for each app, inner ordering is divisions, so dividing by the
577  // number of divisions gets us the index of the application
578  else
579  return kdtree_index / getNumDivisions();
580 }
581 
582 unsigned int
584 {
585  if (_use_nearest_app)
586  return 1;
587  else if (_group_subapps)
588  return _from_meshes.size();
589  else
590  return 1;
591 }
592 
593 unsigned int
595 {
596  // This is not used currently, but conceptually it is better to only divide the domain with the
597  // local of local applications rather than the global number of positions (# global applications
598  // here)
599  if (_use_nearest_app)
600  return _from_meshes.size();
601  // Each nearest-position region is a division
604  EXEC_INITIAL);
605  // Assume all mesh divisions (on each sub-app) has the same number of divisions. This is checked
606  else if (!_from_mesh_divisions.empty())
608  // Grouping subapps or no special mode, we do not subdivide
609  else
610  return 1;
611 }
612 
613 Point
615  const Point & pt) const
616 {
617 
618  if (!_nearest_positions_obj &&
619  (!_from_transforms[getGlobalSourceAppIndex(i_from)]->hasCoordinateSystemTypeChange() ||
621  return _from_transforms[getGlobalSourceAppIndex(i_from)]->mapBack(pt);
623  return pt - _from_positions[i_from];
624  else
625  return pt;
626 }
627 
628 bool
630  const Point & pt, const unsigned int mesh_div, const unsigned int i_from) const
631 {
632  // Only use the KDTree from the closest position if in "nearest-position" mode
634  {
635  // See computeNumSources for the number of sources. i_from is the index in the source loop
636  // i_from is local if looping on _from_problems as sources, positions are indexed globally
637  // i_from is already indexing in positions if using group_subapps
638  auto position_index = i_from; // if _group_subapps
639  if (_use_nearest_app)
640  position_index = getGlobalSourceAppIndex(i_from);
641  else if (!_group_subapps)
642  position_index = i_from % getNumDivisions();
643 
644  // NOTE: if two positions are equi-distant to the point, this will chose one
645  // This problem is detected if using search_value_conflicts in this call
646  if (!closestToPosition(position_index, pt))
647  return false;
648  }
649 
650  // Application index depends on which source/grouping mode we are using
651  const unsigned int app_index = getAppIndex(i_from, i_from / getNumDivisions());
652 
653  // Check mesh restriction before anything
655  {
656  // We have to be careful that getPointInLocalSourceFrame returns in the reference frame
658  mooseError("Nearest-positions + source_app_must_contain_point not implemented");
659  // Transform the point to place it in the local coordinate system
660  const auto local_pt = getPointInLocalSourceFrame(app_index, pt);
661  if (!inMesh(_from_point_locators[app_index].get(), local_pt))
662  return false;
663  }
664 
665  // Check the mesh division. We have handled the restriction of the source locations when
666  // building the nearest-neighbor trees. We only need to check that we meet the required
667  // source division index.
668  if (!_from_mesh_divisions.empty())
669  {
670  mooseAssert(mesh_div != MooseMeshDivision::INVALID_DIVISION_INDEX,
671  "We should not be receiving point requests with an invalid "
672  "source mesh division index");
673  const unsigned int kd_div_index = i_from % getNumDivisions();
674 
675  // If matching source mesh divisions to target apps, we check that the index of the target
676  // application, which was passed in the point request, is equal to the current mesh division
677  if (_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_SUBAPP_INDEX &&
678  mesh_div != kd_div_index)
679  return false;
680  // If matching source mesh divisions to target mesh divisions, we check that the index of the
681  // target mesh division, which was passed in the point request, is equal to the current mesh
682  // division
683  else if ((_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
684  _to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX) &&
685  mesh_div != kd_div_index)
686  return false;
687  }
688 
689  // If matching target apps to source mesh divisions, we check that the global index of the
690  // application is equal to the target mesh division index, which was passed in the point request
691  if (_to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_SUBAPP_INDEX &&
692  mesh_div != getGlobalSourceAppIndex(app_index))
693  return false;
694 
695  return true;
696 }
LAGRANGE
void mooseInfo(Args &&... args) const
Definition: MooseBase.h:317
void renameParam(const std::string &old_name, const std::string &new_name, const std::string &new_docstring)
Rename a parameter and provide a new documentation string.
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)
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.
Definition: MooseUtils.h:380
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
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...
char ** blocks
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.
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
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.
MeshBase & mesh
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
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
FEProblemBase & _fe_problem
Definition: Transfer.h:97
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...
DISCONTINUOUS
auto max(const L &left, const R &right)
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...
An interface that allows the marking of invalid solutions during a solve.
SIDE_DISCONTINUOUS
CONSTANT
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...
Definition: MooseMesh.h:88
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...
Definition: MooseEnum.h:33
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 getNumPositions(bool initial=false) const
}
Definition: Positions.h:37
auto norm(const T &a) -> decltype(std::abs(a))
bool hasSecondOrderElements()
check if the mesh has SECOND order elements
Definition: MooseMesh.C:3711
unsigned int INVALID_DIVISION_INDEX
Invalid subdomain id to return when outside the mesh division.
Definition: MeshDivision.h:28
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.
Definition: MooseBase.h:295
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...
Definition: MooseBase.h:267
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...
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...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:195
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.
uint8_t dof_id_type
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
Definition: Moose.C:30