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