https://mooseframework.inl.gov
MultiAppGeneralFieldTransfer.h
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 
10 #pragma once
11 
13 #include "KDTree.h"
14 #include "PointIndexedMap.h"
15 
16 #include "libmesh/generic_projector.h"
17 #include "libmesh/meshfree_interpolation.h"
18 #include "libmesh/system.h"
19 #include "libmesh/mesh_function.h"
20 #include "libmesh/parallel_algebra.h" // for communicator send and receive stuff
21 
22 class Positions;
23 class MeshDivision;
24 
38 {
39 public:
41 
43 
44  virtual void initialSetup() override;
45  virtual void getAppInfo() override;
46  virtual void execute() override;
47  virtual void postExecute() override;
48 
50  VariableName getFromVarName(unsigned int var_index);
51 
53  VariableName getToVarName(unsigned int var_index);
54 
55 protected:
57  virtual void checkSiblingsTransferSupported() const override {}
58 
59  /*
60  * Prepare evaluation of interpolation values
61  * @param var_index index of the variable & component to prepare for. This routine is called once
62  * for each variable and component to transfer
63  */
64  virtual void prepareEvaluationOfInterpValues(const unsigned int var_index) = 0;
65 
66  /*
67  * Evaluate interpolation values for incoming points
68  * @param incoming_points vector of point requests with an additional integer to add constraints
69  * on the source regions
70  * @param outgoing_vals vector of (evaluated value, distance to value location) for each of the
71  * incoming point requests
72  */
73  virtual void
74  evaluateInterpValues(const std::vector<std::pair<Point, unsigned int>> & incoming_points,
75  std::vector<std::pair<Real, Real>> & outgoing_vals) = 0;
76 
77  /*
78  * Local from bounding boxes for current processor
79  * @param local_bboxes vector of the local (to this process) bounding boxes
80  * for the source applications
81  */
82  void extractLocalFromBoundingBoxes(std::vector<BoundingBox> & local_bboxes);
83 
84  /*
85  * Whether all source mesh checks pass on the given points:
86  * - within the source mesh bounding box
87  * - inside source block restriction
88  * - inside source boundary restriction / in an element near the origin boundary restriction
89  * - inside source mesh division and at the same index as in the target mesh division
90  * - inside app mesh (if not already known to be inside a block or near a boundary)
91  * @param i_from the index of the source problem/mesh
92  * @param local_bboxes the bounding boxes for the local applications
93  * @param pt the point to consider (in reference frame, not source or target frame)
94  * @param mesh_div index of the point to consider in the target mesh division.
95  * If mesh divisions are used to match regions between two meshes, then the index
96  * must match in the target and source mesh divisions
97  * @param distance distance from the point to the target mesh division. This may be used to favor
98  * a point over another, for example in a nearest-positions algorithm where we want
99  * to report the distance of the point pt to the closest position of the 1-NN partition
100  */
101  bool acceptPointInOriginMesh(unsigned int i_from,
102  const std::vector<BoundingBox> & local_bboxes,
103  const Point & pt,
104  const unsigned int mesh_div,
105  Real & distance) const;
106 
107  /*
108  * Whether or not a given point is within the mesh of an origin (from) app
109  * @param pl locator for the mesh of the source app
110  * @param pt point in the local coordinates of the source app we're considering
111  */
112  bool inMesh(const libMesh::PointLocatorBase * const pl, const Point & pt) const;
113 
114  /*
115  * Whether or not a given element is part of the given blocks
116  * Passing the mesh is useful to override the definition of the block restriction
117  */
118  bool inBlocks(const std::set<SubdomainID> & blocks, const Elem * elem) const;
119  virtual bool
120  inBlocks(const std::set<SubdomainID> & blocks, const MooseMesh & mesh, const Elem * elem) const;
121 
122  /*
123  * Whether or not a given node is part of an element in the given blocks
124  * @param blocks list of blocks to examine.
125  * @param mesh containing the element and the node
126  * @param node node to consider
127  * A node is in a block if an element it is part of is in that block
128  */
129  bool
130  inBlocks(const std::set<SubdomainID> & blocks, const MooseMesh & mesh, const Node * node) const;
131 
132  /*
133  * Whether or not a given point is part of an element in the given blocks
134  * @param blocks list of blocks to examine
135  * @param point locator to find the point in the blocks
136  * @param pt point to examine, in the local coordinates (same as the point locator)
137  */
138  bool inBlocks(const std::set<SubdomainID> & blocks,
139  const libMesh::PointLocatorBase * const pl,
140  const Point & pt) const;
141 
142  /*
143  * Whether or not a given node is part of the given boundaries
144  */
145  bool onBoundaries(const std::set<BoundaryID> & boundaries,
146  const MooseMesh & mesh,
147  const Node * node) const;
148 
149  /*
150  * Whether or not a given element is near the specified boundaries
151  * Depending on the '_elemental_boundary_restriction_on_sides' this can mean it shares a side or
152  * a node with the boundary
153  */
154  bool onBoundaries(const std::set<BoundaryID> & boundaries,
155  const MooseMesh & mesh,
156  const Elem * elem) const;
157 
158  /*
159  * Whether or not a point is inside an element that is near the specified boundaries
160  * See onBoundaries(bd, mesh, elem) for definition of near
161  * @param boundaries boundaries of interest for whether the point is near one
162  * @param block_restriction limits the size of the mesh to search for the point.
163  * Note: an empty set means ALL blocks should be considered
164  * @param mesh the mesh to look for the boundaries in
165  * @param pl the point locator that searches the mesh
166  * @param pt the point we want to know whether it is close to a boundary (local coordinates)
167  */
168  bool onBoundaries(const std::set<BoundaryID> & boundaries,
169  const std::set<SubdomainID> & block_restriction,
170  const MooseMesh & mesh,
171  const libMesh::PointLocatorBase * const pl,
172  const Point & pt) const;
173 
182  bool acceptPointMeshDivision(const Point & pt,
183  const unsigned int i_local,
184  const unsigned int only_from_this_mesh_div) const;
185 
192  bool closestToPosition(unsigned int pos_index, const Point & pt) const;
193 
195  const std::vector<unsigned int> _from_var_components;
196 
198  const std::vector<unsigned int> _to_var_components;
199 
203 
205  const bool _use_nearest_app;
206  // NOTE: Keeping track of that distance is not optimal efficiency-wise, because we could find the
207  // closest app and only query the point from there. However, if the value from the closest
208  // app is invalid and the second closest is valid, then the results can vary in parallel
209  // If both apps are the same rank, the closest app is used, the point has an invalid value
210  // If each app are on a different rank, the second closest return a valid value, it gets
211  // used
212 
213  // Positions object to use to match target points and origin points as closest to the same
214  // Position
216 
220 
222  std::set<SubdomainID> _from_blocks;
223 
225  std::set<SubdomainID> _to_blocks;
226 
228  std::set<BoundaryID> _to_boundaries;
229 
231  std::set<BoundaryID> _from_boundaries;
232 
234  std::vector<const MeshDivision *> _from_mesh_divisions;
235 
237  std::vector<const MeshDivision *> _to_mesh_divisions;
238 
241 
244 
247  {
251  };
252 
255 
257  std::vector<std::unique_ptr<libMesh::PointLocatorBase>> _from_point_locators;
258 
261  std::vector<unsigned int> _global_app_start_per_proc;
262 
267 
270 
273 
276 
286  bool detectConflict(Real value_1, Real value_2, Real distance_1, Real distance_2) const;
287 
299  void registerConflict(unsigned int problem, dof_id_type dof_id, Point p, Real dist, bool local);
300 
301 private:
303  std::vector<MooseVariableFieldBase *> _to_variables;
304 
306  typedef std::unordered_map<processor_id_type, std::vector<std::pair<Point, unsigned int>>>
308 
310  struct PointInfo
311  {
312  unsigned int problem_id; // problem id
313  dof_id_type dof_object_id; // node or elem id
314  };
315 
317  struct InterpInfo
318  {
319  processor_id_type pid; // Processor id type
320  Real interp; // Interpolation
321  Real distance; // distance from target to source
322  };
323 
325  typedef std::unordered_map<processor_id_type, std::vector<PointInfo>> ProcessorToPointInfoVec;
326 
328  typedef std::vector<std::unordered_map<dof_id_type, InterpInfo>> DofobjectToInterpValVec;
329 
334 
336  typedef std::vector<InterpCache> InterpCaches;
337 
339  unsigned int _var_size;
340 
343 
346 
349 
351  std::vector<Real> _fixed_bbox_size;
352 
354  std::vector<unsigned int> _froms_per_proc;
355 
360  std::vector<BoundingBox> _from_bboxes;
361 
364 
369  std::vector<std::tuple<unsigned int, dof_id_type, Point, Real>> _local_conflicts;
370 
375  std::vector<std::tuple<unsigned int, dof_id_type, Point, Real>> _received_conflicts;
376 
380  void prepareToTransfer();
381 
385  void transferVariable(unsigned int i);
386 
387  /*
388  * Extract to-points for which we need to compute interpolation values on the source domains
389  */
390  void extractOutgoingPoints(const unsigned int var_index, ProcessorToPointVec & outgoing_points);
391 
392  /*
393  * Which processors include this point
394  */
395  void locatePointReceivers(const Point point, std::set<processor_id_type> & processors);
396 
397  /*
398  * cache incoming values
399  */
401  processor_id_type pid,
402  const unsigned int var_index,
403  std::vector<PointInfo> & pointInfoVec,
404  const std::vector<std::pair<Point, unsigned int>> & point_requests,
405  const std::vector<std::pair<Real, Real>> & incoming_vals,
406  DofobjectToInterpValVec & dofobject_to_valsvec, // for nodal + constant monomial
407  InterpCaches & interp_caches, // for higher order elemental values
408  InterpCaches & distance_caches); // same but helps make origin point decisions
409 
419  void examineReceivedValueConflicts(const unsigned int var_index,
420  const DofobjectToInterpValVec & dofobject_to_valsvec,
421  const InterpCaches & distance_caches);
422 
433  void examineLocalValueConflicts(const unsigned int var_index,
434  const DofobjectToInterpValVec & dofobject_to_valsvec,
435  const InterpCaches & distance_caches);
436 
438  void outputValueConflicts(const unsigned int var_index,
439  const DofobjectToInterpValVec & dofobject_to_valsvec,
440  const InterpCaches & distance_caches);
441 
442  /*
443  * Set values to solution
444  * @param var the variable to set
445  * @param dofobject_to_valsvec a vector of maps from DoF to values, for each to_problem
446  * Used for nodal + constant monomial variables
447  * @param interp_caches a vector of maps from point to value, for each to_problem
448  * Used for higher order elemental variables
449  */
450  void setSolutionVectorValues(const unsigned int var_index,
451  const DofobjectToInterpValVec & dofobject_to_valsvec,
452  const InterpCaches & interp_caches);
453 
454  /*
455  * Cache pointInfo
456  */
457  void cacheOutgoingPointInfo(const Point point,
458  const dof_id_type dof_object_id,
459  const unsigned int problem_id,
460  ProcessorToPointVec & outgoing_points);
461 
467  Real bboxMinDistance(const Point & p, const BoundingBox & bbox) const;
468 
474  Real bboxMaxDistance(const Point & p, const BoundingBox & bbox) const;
475 
480  Point getMaxToProblemsBBoxDimensions() const;
481 
485  std::vector<BoundingBox> getRestrictedFromBoundingBoxes() const;
486 
491  std::vector<unsigned int> getGlobalStartAppPerProc() const;
492 };
493 
494 // Anonymous namespace for data, functors to use with GenericProjector.
496 {
497 // Transfer::OutOfMeshValue is an actual number. Why? Why!
498 static_assert(std::numeric_limits<Real>::has_infinity,
499  "What are you trying to use for Real? It lacks infinity!");
501 
502 inline bool
504 {
505  // Might need to be changed for e.g. NaN
507 }
508 
509 // We need two functors that record point (value and gradient,
510 // respectively) requests, so we know what queries we need to make
511 // to other processors
512 
516 template <typename Output>
518 {
519 protected:
521 
522 public:
526 
528 
529  RecordRequests(RecordRequests & primary) : _primary(&primary) {}
530 
532  {
533  if (_primary)
534  {
535  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
536  _primary->_points_requested.insert(
538  }
539  }
540 
542 
544  unsigned int /*variable_index*/,
545  unsigned int /*elem_dim*/,
546  const Node & n,
547  bool /*extra_hanging_dofs*/,
548  const Real /*time*/)
549  {
550  _points_requested.push_back(n);
551  return 0;
552  }
553 
555  unsigned int /*variable_index*/,
556  const Point & n,
557  const Real /*time*/,
558  bool /*skip_context_check*/)
559  {
560  _points_requested.push_back(n);
561  return 0;
562  }
563 
564  bool is_grid_projection() { return false; }
565 
567  unsigned int /*i*/,
568  unsigned int /*dim*/,
569  const Node & /*n*/,
570  std::vector<Output> & /*derivs*/)
571  {
572  mooseError("Not implemented");
573  } // this is only for grid projections
574 
576  const Elem &, unsigned int, unsigned int, std::vector<dof_id_type> &, std::vector<Output> &)
577  {
578  mooseError("Not implemented");
579  }
580 
581  void eval_old_dofs(const Elem &,
582  const libMesh::FEType &,
583  unsigned int,
584  unsigned int,
585  std::vector<dof_id_type> &,
586  std::vector<Output> &)
587  {
588  mooseError("Not implemented");
589  }
590 
591  std::vector<Point> & points_requested() { return _points_requested; }
592 
593 private:
595  std::vector<Point> _points_requested;
596 
598 };
599 
600 // We need a null action functor to use
601 // with them (because we won't be ready to set any values at that point)
602 template <typename Val>
604 {
605 public:
606  typedef Val InsertInput;
607 
609 
610  void insert(dof_id_type, Val) {}
611 
612  void insert(const std::vector<dof_id_type> &, const DenseVector<Val> &) {}
613 };
614 
615 // We need two functors that respond to point (value and gradient,
616 // respectively) requests based on the cached values of queries answered by
617 // other processors.
618 
622 template <typename Output>
624 {
625 protected:
627 
628 public:
630 
634 
640  CachedData(const Cache & cache, const libMesh::FunctionBase<Output> & backup, Real default_value)
641  : _cache(cache), _backup(backup.clone()), _default_value(default_value)
642  {
643  }
644 
646  CachedData(const CachedData & primary)
647  : _cache(primary._cache),
648  _backup(primary._backup->clone()),
650  {
651  }
652 
654 
657  unsigned int /*i*/,
658  unsigned int /*elem_dim*/,
659  const Node & n,
660  bool /*extra_hanging_dofs*/,
661  const Real /*time*/)
662  {
663  auto it = _cache.find(n);
664  if (it == _cache.end())
665  {
667  return _default_value;
668  else
669  return (*_backup)(n);
670  }
671  else
672  return it->second;
673  }
674 
677  unsigned int /*i*/,
678  const Point & n,
679  const Real /*time*/,
680  bool /*skip_context_check*/)
681  {
682  auto it = _cache.find(n);
683  if (it == _cache.end())
684  {
686  return _default_value;
687  else
688  return (*_backup)(n);
689  }
690  else
691  return it->second;
692  }
693 
694  bool is_grid_projection() { return false; }
695 
697  unsigned int /*i*/,
698  unsigned int /*dim*/,
699  const Node & /*n*/,
700  std::vector<Output> & /*derivs*/)
701  {
702  mooseError("Not implemented");
703  } // this is only for grid projections
704 
706  const Elem &, unsigned int, unsigned int, std::vector<dof_id_type> &, std::vector<Output> &)
707  {
708  mooseError("Not implemented");
709  }
710 
711  void eval_old_dofs(const Elem &,
712  const libMesh::FEType &,
713  unsigned int,
714  unsigned int,
715  std::vector<dof_id_type> &,
716  std::vector<Output> &)
717  {
718  mooseError("Not implemented");
719  }
720 
721 private:
723  const Cache & _cache;
724 
726  std::unique_ptr<libMesh::FunctionBase<Output>> _backup;
727 
730 };
731 
732 }
std::vector< BoundingBox > _from_bboxes
Bounding boxes for all source applications.
virtual void evaluateInterpValues(const std::vector< std::pair< Point, unsigned int >> &incoming_points, std::vector< std::pair< Real, Real >> &outgoing_vals)=0
void examineReceivedValueConflicts(const unsigned int var_index, const DofobjectToInterpValVec &dofobject_to_valsvec, const InterpCaches &distance_caches)
Remove potential value conflicts that did not materialize because another source was closer Several e...
void cacheOutgoingPointInfo(const Point point, const dof_id_type dof_object_id, const unsigned int problem_id, ProcessorToPointVec &outgoing_points)
bool _source_app_must_contain_point
Whether the source app mesh must actually contain the points for them to be considered or whether the...
char ** blocks
void registerConflict(unsigned int problem, dof_id_type dof_id, Point p, Real dist, bool local)
Register a potential value conflict, e.g.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
std::unordered_map< Point, Number, hash_point >::const_iterator end() const
void eval_mixed_derivatives(const libMesh::FEMContext &, unsigned int, unsigned int, const Node &, std::vector< Output > &)
std::vector< Point > _points_requested
Vector of points requested.
VariableName getFromVarName(unsigned int var_index)
Get the source variable name, with the suffix for array/vector variables.
libMesh::TensorTools::MakeReal< Output >::type RealType
Positions objects are under the hood Reporters.
Definition: Positions.h:20
void locatePointReceivers(const Point point, std::set< processor_id_type > &processors)
std::vector< std::unique_ptr< libMesh::PointLocatorBase > > _from_point_locators
Point locators, useful to examine point location with regards to domain restriction.
virtual void execute() override
Execute the transfer.
MultiAppGeneralFieldTransfer(const InputParameters &parameters)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
bool acceptPointMeshDivision(const Point &pt, const unsigned int i_local, const unsigned int only_from_this_mesh_div) const
Whether a point lies inside the mesh division delineated by the MeshDivision object.
std::vector< unsigned int > getGlobalStartAppPerProc() const
Get global index for the first app each processes owns Requires a global communication, must be called on every domain simultaneously.
Output eval_at_point(const libMesh::FEMContext &, unsigned int, const Point &n, const Real, bool)
Gets a value at a point.
std::vector< Real > _fixed_bbox_size
Set the bounding box sizes manually.
std::vector< const MeshDivision * > _from_mesh_divisions
Division of the origin mesh.
std::vector< const MeshDivision * > _to_mesh_divisions
Division of the target mesh.
unsigned int _var_size
The number of variables to transfer.
Base class for MeshDivision objects.
Definition: MeshDivision.h:35
const Real _default_extrapolation_value
Value to use when no received data is valid for a target location.
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.
const std::vector< unsigned int > _from_var_components
Origin array/vector variable components.
void eval_old_dofs(const Elem &, unsigned int, unsigned int, std::vector< dof_id_type > &, std::vector< Output > &)
virtual void getAppInfo() override
This method will fill information into the convenience member variables (_to_problems, _from_meshes, etc.)
VariableName getToVarName(unsigned int var_index)
Get the target variable name, with the suffix for array/vector variables.
std::set< BoundaryID > _from_boundaries
Origin boundary(ies) restriction.
void outputValueConflicts(const unsigned int var_index, const DofobjectToInterpValVec &dofobject_to_valsvec, const InterpCaches &distance_caches)
Report on conflicts between overlapping child apps, equidistant origin points etc.
std::unordered_map< processor_id_type, std::vector< PointInfo > > ProcessorToPointInfoVec
A map from pid to a set of point info.
Real bboxMinDistance(const Point &p, const BoundingBox &bbox) const
Compute minimum distance.
PointIndexedMap InterpCache
A map from Point to interpolation values NOTE: this is not an asynchronous cache. ...
void extractOutgoingPoints(const unsigned int var_index, ProcessorToPointVec &outgoing_points)
std::set< SubdomainID > _to_blocks
Target block(s) restriction.
Based class for output objects.
Definition: Output.h:43
std::unordered_map< processor_id_type, std::vector< std::pair< Point, unsigned int > > > ProcessorToPointVec
A map from pid to a set of points.
Point getMaxToProblemsBBoxDimensions() const
Obtains the max dimensions to scale all points in the mesh.
const Cache & _cache
Data to return for cached points.
uint8_t processor_id_type
virtual void prepareEvaluationOfInterpValues(const unsigned int var_index)=0
const std::vector< unsigned int > _to_var_components
Target array/vector variable components.
Real bboxMaxDistance(const Point &p, const BoundingBox &bbox) const
Compute max distance.
std::vector< std::tuple< unsigned int, dof_id_type, Point, Real > > _local_conflicts
Keeps track of all local equidistant points to requested points, creating an indetermination in which...
const bool _use_bounding_boxes
Whether to use bounding boxes to determine the applications that may receive point requests then send...
libMesh::TensorTools::MakeBaseNumber< Output >::type DofValueType
Value request recording base class.
bool acceptPointInOriginMesh(unsigned int i_from, const std::vector< BoundingBox > &local_bboxes, const Point &pt, const unsigned int mesh_div, Real &distance) const
CachedData(const Cache &cache, const libMesh::FunctionBase< Output > &backup, Real default_value)
Constructor.
std::set< BoundaryID > _to_boundaries
Target boundary(ies) restriction.
std::vector< unsigned int > _froms_per_proc
Number of source/from applications per processor. This vector is indexed by processor id...
std::vector< std::tuple< unsigned int, dof_id_type, Point, Real > > _received_conflicts
Keeps track of all received conflicts.
MeshDivisionTransferUse
Matching enum for the mesh division behaviors.
std::vector< std::unordered_map< dof_id_type, InterpInfo > > DofobjectToInterpValVec
A vector, indexed by to-problem id, of maps from dof object to interpolation values.
Real _bbox_factor
How much we should relax bounding boxes.
bool _error_on_miss
Error out when some points can not be located.
libMesh::TensorTools::MakeBaseNumber< Output >::type DofValueType
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:88
virtual void checkSiblingsTransferSupported() const override
Siblings transfers fully supported.
Transfers variables on possibly different meshes while conserving a user defined property (Postproces...
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::unique_ptr< libMesh::FunctionBase< Output > > _backup
Function to evaluate for uncached points.
std::unordered_map< Point, Number, hash_point >::const_iterator find(const Point &pt) const
An unordered map indexed by Point, eg 3 floating point numbers Because floating point rounding errors...
void transferVariable(unsigned int i)
Performs the transfer for the variable of index i.
Output eval_at_node(const libMesh::FEMContext &, unsigned int, unsigned int, const Node &n, bool, const Real)
Gets a value at the node location.
void examineLocalValueConflicts(const unsigned int var_index, const DofobjectToInterpValVec &dofobject_to_valsvec, const InterpCaches &distance_caches)
Remove potential value conflicts that did not materialize because another source was closer Several e...
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.
std::vector< BoundingBox > getRestrictedFromBoundingBoxes() const
Get from bounding boxes for given domains and boundaries.
libMesh::TensorTools::MakeReal< Output >::type RealType
void eval_old_dofs(const Elem &, const libMesh::FEType &, unsigned int, unsigned int, std::vector< dof_id_type > &, std::vector< Output > &)
void insert(const std::vector< dof_id_type > &, const DenseVector< Val > &)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool _already_output_search_value_conflicts
Whether we already output the search value conflicts.
const MooseEnum & _to_mesh_division_behavior
How to use the target mesh divisions to restrict the transfer.
virtual void postExecute() override
Add some extra work if necessary after execute().
const Real _default_value
Default value when no point is found.
Value request response base class.
std::set< SubdomainID > _from_blocks
Origin block(s) restriction.
const unsigned int _search_value_conflicts_max_log
How many conflicts are output to console.
std::vector< InterpCache > InterpCaches
A vector of such caches, indexed by to_problem.
void extractLocalFromBoundingBoxes(std::vector< BoundingBox > &local_bboxes)
const InputParameters & parameters() const
Get the parameters of the object.
void eval_old_dofs(const Elem &, unsigned int, unsigned int, std::vector< dof_id_type > &, std::vector< Output > &)
void eval_old_dofs(const Elem &, const libMesh::FEType &, unsigned int, unsigned int, std::vector< dof_id_type > &, std::vector< Output > &)
std::vector< MooseVariableFieldBase * > _to_variables
The target variables.
Real Number
ProcessorToPointInfoVec _processor_to_pointInfoVec
A map from processor to pointInfo vector.
const bool _elemental_boundary_restriction_on_sides
Whether elemental variable boundary restriction is considered by element side or element nodes...
virtual void initialSetup() override
Method called at the beginning of the simulation for checking integrity or doing one-time setup...
void prepareToTransfer()
Initialize supporting attributes like bounding boxes, processor app indexes etc.
CachedData(const CachedData &primary)
Copy constructor.
void eval_mixed_derivatives(const libMesh::FEMContext &, unsigned int, unsigned int, const Node &, std::vector< Output > &)
Output eval_at_node(const libMesh::FEMContext &, unsigned int, unsigned int, const Node &n, bool, const Real)
void setSolutionVectorValues(const unsigned int var_index, const DofobjectToInterpValVec &dofobject_to_valsvec, const InterpCaches &interp_caches)
bool onBoundaries(const std::set< BoundaryID > &boundaries, const MooseMesh &mesh, const Node *node) const
Output eval_at_point(const libMesh::FEMContext &, unsigned int, const Point &n, const Real, bool)
It is a general field transfer.
bool _greedy_search
Whether or not a greedy strategy will be used If true, all the partitions will be checked for a given...
void cacheIncomingInterpVals(processor_id_type pid, const unsigned int var_index, std::vector< PointInfo > &pointInfoVec, const std::vector< std::pair< Point, unsigned int >> &point_requests, const std::vector< std::pair< Real, Real >> &incoming_vals, DofobjectToInterpValVec &dofobject_to_valsvec, InterpCaches &interp_caches, InterpCaches &distance_caches)
std::vector< unsigned int > _global_app_start_per_proc
First app each processor owns, indexed by processor If no app on the processor, will have a -1 for th...
uint8_t dof_id_type
const MooseEnum & _from_mesh_division_behavior
How to use the origin mesh divisions to restrict the transfer.
bool inBlocks(const std::set< SubdomainID > &blocks, const Elem *elem) const
const bool _use_nearest_app
Whether to keep track of the distance from the requested point to the app position.