Line data Source code
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 :
12 : #include "MultiAppConservativeTransfer.h"
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 :
25 : /**
26 : * It is a general field transfer. It will do the following things
27 : * 1) From part of source domain to part of domain. Support subdomains/boundaries to
28 : * subdomains/boundaries, mixing as appropriate
29 : * 2) interpolation and extrapolation, as appropriate
30 : * 3) Support higher order FEM
31 : * 4) Support mixed orders between source and target variables
32 : * 5) Support both distributed and replicated meshes
33 : * 6) Support both origin and target displaced meshes
34 : * 7) Support siblings transfers
35 : * 8) Support multiple child apps in both the transfer source and target
36 : */
37 : class MultiAppGeneralFieldTransfer : public MultiAppConservativeTransfer
38 : {
39 : public:
40 : static InputParameters validParams();
41 :
42 : MultiAppGeneralFieldTransfer(const InputParameters & parameters);
43 :
44 : virtual void initialSetup() override;
45 : virtual void getAppInfo() override;
46 : virtual void execute() override;
47 : virtual void postExecute() override;
48 :
49 : /// Get the source variable name, with the suffix for array/vector variables
50 : VariableName getFromVarName(unsigned int var_index) const;
51 :
52 : /// Get the target variable name, with the suffix for array/vector variables
53 : VariableName getToVarName(unsigned int var_index);
54 :
55 : protected:
56 : /// Siblings transfers fully supported
57 1276 : virtual void checkSiblingsTransferSupported() const override {}
58 :
59 : /// Return a pointer to a target variable
60 109576 : MooseVariableFieldBase * getToVariable(unsigned int var_index) const
61 : {
62 109576 : return _to_variables[var_index];
63 : }
64 :
65 : /**
66 : * Return a human-readable description of the data source (variable, functor, user object, etc.)
67 : * used for conflict warning messages. Override in derived classes that use a different source
68 : * type (e.g. functors).
69 : * @param var_index index of the variable/functor being transferred
70 : */
71 : virtual std::string getDataSourceName(unsigned int var_index) const;
72 :
73 : /*
74 : * Prepare evaluation of interpolation values
75 : * @param var_index index of the variable & component to prepare for. This routine is called once
76 : * for each variable and component to transfer
77 : */
78 : virtual void prepareEvaluationOfInterpValues(const unsigned int var_index) = 0;
79 :
80 : /*
81 : * Evaluate interpolation values for incoming points
82 : * @param var_index the index of the variable being transferred (same for source & target)
83 : * @param incoming_points vector of point requests with an additional integer to add constraints
84 : * on the source regions
85 : * @param outgoing_vals vector of (evaluated value, distance to value location) for each of the
86 : * incoming point requests
87 : */
88 : virtual void
89 : evaluateInterpValues(const unsigned int var_index,
90 : const std::vector<std::pair<Point, unsigned int>> & incoming_points,
91 : std::vector<std::pair<Real, Real>> & outgoing_vals) = 0;
92 :
93 : /*
94 : * Local from bounding boxes for current processor
95 : * @param local_bboxes vector of the local (to this process) bounding boxes
96 : * for the source applications
97 : */
98 : void extractLocalFromBoundingBoxes(std::vector<BoundingBox> & local_bboxes);
99 :
100 : /*
101 : * Whether all source mesh checks pass on the given points:
102 : * - within the source mesh bounding box
103 : * - inside source block restriction
104 : * - inside source boundary restriction / in an element near the origin boundary restriction
105 : * - inside source mesh division and at the same index as in the target mesh division
106 : * - inside app mesh (if not already known to be inside a block or near a boundary)
107 : * @param i_from the index of the source problem/mesh
108 : * @param local_bboxes the bounding boxes for the local applications
109 : * @param pt the point to consider (in reference frame, not source or target frame)
110 : * @param mesh_div index of the point to consider in the target mesh division.
111 : * If mesh divisions are used to match regions between two meshes, then the index
112 : * must match in the target and source mesh divisions
113 : * @param distance distance from the point to the target mesh division. This may be used to favor
114 : * a point over another, for example in a nearest-positions algorithm where we want
115 : * to report the distance of the point pt to the closest position of the 1-NN partition
116 : */
117 : bool acceptPointInOriginMesh(unsigned int i_from,
118 : const std::vector<BoundingBox> & local_bboxes,
119 : const Point & pt,
120 : const unsigned int mesh_div,
121 : Real & distance) const;
122 :
123 : /*
124 : * Whether or not a given point is within the mesh of an origin (from) app
125 : * @param pl locator for the mesh of the source app
126 : * @param pt point in the local coordinates of the source app we're considering
127 : */
128 : bool inMesh(const libMesh::PointLocatorBase * const pl, const Point & pt) const;
129 :
130 : /*
131 : * Whether or not a given element is part of the given blocks
132 : * Passing the mesh is useful to override the definition of the block restriction
133 : */
134 : bool inBlocks(const std::set<SubdomainID> & blocks, const Elem * elem) const;
135 : virtual bool
136 : inBlocks(const std::set<SubdomainID> & blocks, const MooseMesh & mesh, const Elem * elem) const;
137 :
138 : /*
139 : * Whether or not a given node is part of an element in the given blocks
140 : * @param blocks list of blocks to examine.
141 : * @param mesh containing the element and the node
142 : * @param node node to consider
143 : * A node is in a block if an element it is part of is in that block
144 : */
145 : bool
146 : inBlocks(const std::set<SubdomainID> & blocks, const MooseMesh & mesh, const Node * node) const;
147 :
148 : /*
149 : * Whether or not a given point is part of an element in the given blocks
150 : * @param blocks list of blocks to examine
151 : * @param point locator to find the point in the blocks
152 : * @param pt point to examine, in the local coordinates (same as the point locator)
153 : */
154 : bool inBlocks(const std::set<SubdomainID> & blocks,
155 : const libMesh::PointLocatorBase * const pl,
156 : const Point & pt) const;
157 :
158 : /*
159 : * Whether or not a given node is part of the given boundaries
160 : */
161 : bool onBoundaries(const std::set<BoundaryID> & boundaries,
162 : const MooseMesh & mesh,
163 : const Node * node) const;
164 :
165 : /*
166 : * Whether or not a given element is near the specified boundaries
167 : * Depending on the '_elemental_boundary_restriction_on_sides' this can mean it shares a side or
168 : * a node with the boundary
169 : */
170 : bool onBoundaries(const std::set<BoundaryID> & boundaries,
171 : const MooseMesh & mesh,
172 : const Elem * elem) const;
173 :
174 : /*
175 : * Whether or not a point is inside an element that is near the specified boundaries
176 : * See onBoundaries(bd, mesh, elem) for definition of near
177 : * @param boundaries boundaries of interest for whether the point is near one
178 : * @param block_restriction limits the size of the mesh to search for the point.
179 : * Note: an empty set means ALL blocks should be considered
180 : * @param mesh the mesh to look for the boundaries in
181 : * @param pl the point locator that searches the mesh
182 : * @param pt the point we want to know whether it is close to a boundary (local coordinates)
183 : */
184 : bool onBoundaries(const std::set<BoundaryID> & boundaries,
185 : const std::set<SubdomainID> & block_restriction,
186 : const MooseMesh & mesh,
187 : const libMesh::PointLocatorBase * const pl,
188 : const Point & pt) const;
189 :
190 : /**
191 : * Whether a point lies inside the mesh division delineated by the MeshDivision object
192 : * @param pt point to examine, in the local coordinates (source frame for from_direction=true)
193 : * @param i_local the index of the problem to consider, holding the mesh division to examine
194 : * @param only_from_this_mesh_div a mesh division index that must be matched when the
195 : * to/from_mesh_division_behavior for the direction examined is MATCH_DIVISION/SUBAPP_INDEX
196 : * It is ignored otherwise
197 : */
198 : bool acceptPointMeshDivision(const Point & pt,
199 : const unsigned int i_local,
200 : const unsigned int only_from_this_mesh_div) const;
201 :
202 : /**
203 : * Whether a point is closest to a position at the index specified than any other position
204 : * @param pos_index the index of the position to consider in the positions vector
205 : * @param pt the point
206 : * @return whether the point is closest to this position than any other in the positions vector
207 : */
208 : bool closestToPosition(unsigned int pos_index, const Point & pt) const;
209 :
210 : /// Origin array/vector variable components
211 : const std::vector<unsigned int> _from_var_components;
212 :
213 : /// Target array/vector variable components
214 : const std::vector<unsigned int> _to_var_components;
215 :
216 : /// Whether to use bounding boxes to determine the applications that may receive point requests
217 : /// then send value data, and at other various checks
218 : const bool _use_bounding_boxes;
219 :
220 : /// Whether to keep track of the distance from the requested point to the app position
221 : const bool _use_nearest_app;
222 : // NOTE: Keeping track of that distance is not optimal efficiency-wise, because we could find the
223 : // closest app and only query the point from there. However, if the value from the closest
224 : // app is invalid and the second closest is valid, then the results can vary in parallel
225 : // If both apps are the same rank, the closest app is used, the point has an invalid value
226 : // If each app are on a different rank, the second closest return a valid value, it gets
227 : // used
228 :
229 : // Positions object to use to match target points and origin points as closest to the same
230 : // Position
231 : const Positions * _nearest_positions_obj;
232 :
233 : /// Whether the source app mesh must actually contain the points for them to be considered or whether
234 : /// the bounding box is enough. If false, we can interpolate between apps
235 : bool _source_app_must_contain_point;
236 :
237 : /// Origin block(s) restriction
238 : std::set<SubdomainID> _from_blocks;
239 :
240 : /// Target block(s) restriction
241 : std::set<SubdomainID> _to_blocks;
242 :
243 : /// Target boundary(ies) restriction
244 : std::set<BoundaryID> _to_boundaries;
245 :
246 : /// Origin boundary(ies) restriction
247 : std::set<BoundaryID> _from_boundaries;
248 :
249 : /// Division of the origin mesh
250 : std::vector<const MeshDivision *> _from_mesh_divisions;
251 :
252 : /// Division of the target mesh
253 : std::vector<const MeshDivision *> _to_mesh_divisions;
254 :
255 : /// How to use the origin mesh divisions to restrict the transfer
256 : const MooseEnum & _from_mesh_division_behavior;
257 :
258 : /// How to use the target mesh divisions to restrict the transfer
259 : const MooseEnum & _to_mesh_division_behavior;
260 :
261 : /// Matching enum for the mesh division behaviors.
262 : enum MeshDivisionTransferUse
263 : {
264 : RESTRICTION,
265 : MATCH_DIVISION_INDEX,
266 : MATCH_SUBAPP_INDEX
267 : };
268 :
269 : /// Whether elemental variable boundary restriction is considered by element side or element nodes
270 : const bool _elemental_boundary_restriction_on_sides;
271 :
272 : /// Point locators, useful to examine point location with regards to domain restriction
273 : std::vector<std::unique_ptr<libMesh::PointLocatorBase>> _from_point_locators;
274 :
275 : /// First app each processor owns, indexed by processor
276 : /// If no app on the processor, will have a -1 for the app start instead
277 : std::vector<unsigned int> _global_app_start_per_proc;
278 :
279 : /// Whether or not a greedy strategy will be used
280 : /// If true, all the partitions will be checked for a given
281 : /// outgoing point
282 : bool _greedy_search;
283 :
284 : /// Whether to look for conflicts between origin points, multiple valid values for a target point
285 : bool _search_value_conflicts;
286 :
287 : /// Whether we already output the search value conflicts
288 : bool _already_output_search_value_conflicts;
289 :
290 : /// How many conflicts are output to console
291 : const unsigned int _search_value_conflicts_max_log;
292 :
293 : /// How to post treat after the transfer
294 : const MooseEnum _post_transfer_extrapolation;
295 :
296 : /**
297 : * @brief Detects whether two source values are valid and equidistant for a desired target
298 : * location
299 : * @param value_1 value from the first value source / subapp
300 : * @param value_2 value from the second value source / subapp
301 : * @param distance_1 distance from the first source
302 : * @param distance_2 distance from the second source
303 : * @return true if the values are different and distances from the source points/apps are the same
304 : */
305 : bool detectConflict(Real value_1, Real value_2, Real distance_1, Real distance_2) const;
306 :
307 : /**
308 : * Register a potential value conflict, e.g. two or more equidistant source points for a single
309 : * target point, with different values possible
310 : * @param problem problem ID for the point of interest.
311 : * For local conflicts, use origin problem id, for received conflicts, use target id
312 : * @param dof_id id id of the DoF is transferring a DoF. If not, use -1
313 : * @param p point where the conflict happens
314 : * @param dist distance between the origin and the target
315 : * @param local if true, local conflict found when gathering data to send, if false,
316 : * received value conflict found when receiving data from multiple source problems
317 : */
318 : void registerConflict(unsigned int problem, dof_id_type dof_id, Point p, Real dist, bool local);
319 :
320 : private:
321 : /// The target variables
322 : std::vector<MooseVariableFieldBase *> _to_variables;
323 :
324 : /// A map from pid to a set of points
325 : typedef std::unordered_map<processor_id_type, std::vector<std::pair<Point, unsigned int>>>
326 : ProcessorToPointVec;
327 :
328 : /// Point information
329 : struct PointInfo
330 : {
331 : unsigned int problem_id; // problem id
332 : dof_id_type dof_object_id; // node or elem id
333 : };
334 :
335 : /// InterpInfo
336 : struct InterpInfo
337 : {
338 : processor_id_type pid; // Processor id type
339 : Real interp; // Interpolation
340 : Real distance; // distance from target to source
341 : };
342 :
343 : /// A map from pid to a set of point info
344 : typedef std::unordered_map<processor_id_type, std::vector<PointInfo>> ProcessorToPointInfoVec;
345 :
346 : /// A vector, indexed by to-problem id, of maps from dof object to interpolation values
347 : typedef std::vector<std::unordered_map<dof_id_type, InterpInfo>> DofobjectToInterpValVec;
348 :
349 : /// A map from Point to interpolation values
350 : /// NOTE: this is not an asynchronous cache. It is built to completion during the transfer
351 : /// and used as a whole to reconstruct the target variable
352 : typedef PointIndexedMap InterpCache;
353 :
354 : /// A vector of such caches, indexed by to_problem
355 : typedef std::vector<InterpCache> InterpCaches;
356 :
357 : /// The number of variables to transfer
358 : unsigned int _var_size;
359 :
360 : /// Error out when some points can not be located
361 : bool _error_on_miss;
362 :
363 : /// Value to use when no received data is valid for a target location
364 : const Real _default_extrapolation_value;
365 :
366 : /// How much we should relax bounding boxes
367 : Real _bbox_factor;
368 :
369 : /// Set the bounding box sizes manually
370 : std::vector<Real> _fixed_bbox_size;
371 :
372 : /// Number of source/from applications per processor. This vector is indexed by processor id
373 : std::vector<unsigned int> _froms_per_proc;
374 :
375 : /// Bounding boxes for all source applications. The indexing of this vector is similar to
376 : ///'processor_id * n_local_subapps + i_local_subapp', except the number of local subapps can
377 : /// be different on each processor. Use the _from_per_proc vector to find the start/end index for a given
378 : /// processor. See MultiAppGeneralFieldTransfer::locatePointReceivers() for an example
379 : std::vector<BoundingBox> _from_bboxes;
380 :
381 : /// A map from processor to pointInfo vector
382 : ProcessorToPointInfoVec _processor_to_pointInfoVec;
383 :
384 : /// Keeps track of all local equidistant points to requested points, creating an indetermination
385 : /// in which values should be sent for that request
386 : /// We keep the origin problem ID, the dof ID, the point, and the distance origin-target
387 : /// If using nearest-positions the origin problem ID is not set
388 : std::vector<std::tuple<unsigned int, dof_id_type, Point, Real>> _local_conflicts;
389 :
390 : /// Keeps track of all received conflicts. Multiple problems (different subapps for example)
391 : /// are sending values for a target point that do not match and are equally valid/distant
392 : /// We keep the target problem ID, the point/dof ID, the point, and the origin-target distance.
393 : /// The distance indicates whether a potential conflict ended up materializing
394 : std::vector<std::tuple<unsigned int, dof_id_type, Point, Real>> _received_conflicts;
395 :
396 : /**
397 : * Initialize supporting attributes like bounding boxes, processor app indexes etc
398 : */
399 : void prepareToTransfer();
400 :
401 : /**
402 : * Performs the transfer for the variable of index i
403 : */
404 : void transferVariable(unsigned int i);
405 :
406 : /*
407 : * Extract to-points for which we need to compute interpolation values on the source domains
408 : */
409 : void extractOutgoingPoints(const unsigned int var_index, ProcessorToPointVec & outgoing_points);
410 :
411 : /*
412 : * Which processors include this point
413 : */
414 : void locatePointReceivers(const Point point, std::set<processor_id_type> & processors);
415 :
416 : /*
417 : * cache incoming values
418 : */
419 : void cacheIncomingInterpVals(
420 : processor_id_type pid,
421 : const unsigned int var_index,
422 : std::vector<PointInfo> & pointInfoVec,
423 : const std::vector<std::pair<Point, unsigned int>> & point_requests,
424 : const std::vector<std::pair<Real, Real>> & incoming_vals,
425 : DofobjectToInterpValVec & dofobject_to_valsvec, // for nodal + constant monomial
426 : InterpCaches & interp_caches, // for higher order elemental values
427 : InterpCaches & distance_caches); // same but helps make origin point decisions
428 :
429 : /**
430 : * Remove potential value conflicts that did not materialize because another source was closer
431 : * Several equidistant valid values were received, but they were not closest
432 : * @param var_index the index of the variable of interest
433 : * @param dofobject_to_valsvec a data structure mapping dofobjects to received values and
434 : * distances (used for nodal-value-dof-only variables and constant monomials)
435 : * @param distance_caches a cache holding the distances received (used for higher order elemental
436 : * variables)
437 : */
438 : void examineReceivedValueConflicts(const unsigned int var_index,
439 : const DofobjectToInterpValVec & dofobject_to_valsvec,
440 : const InterpCaches & distance_caches);
441 :
442 : /**
443 : * Remove potential value conflicts that did not materialize because another source was closer
444 : * Several equidistant valid values were found when computing values to send, but they were not
445 : * closest, another value got selected
446 : * @param var_index the index of the variable of interest
447 : * @param dofobject_to_valsvec a data structure mapping dofobjects to received values and
448 : * distances (used for nodal-value-dof-only variables and constant monomials)
449 : * @param distance_caches a cache holding the distances received (used for higher order elemental
450 : * variables)
451 : */
452 : void examineLocalValueConflicts(const unsigned int var_index,
453 : const DofobjectToInterpValVec & dofobject_to_valsvec,
454 : const InterpCaches & distance_caches);
455 :
456 : /// Report on conflicts between overlapping child apps, equidistant origin points etc
457 : void outputValueConflicts(const unsigned int var_index,
458 : const DofobjectToInterpValVec & dofobject_to_valsvec,
459 : const InterpCaches & distance_caches);
460 :
461 : /*
462 : * Set values to solution
463 : * @param var the variable to set
464 : * @param dofobject_to_valsvec a vector of maps from DoF to values, for each to_problem
465 : * Used for nodal + constant monomial variables
466 : * @param interp_caches a vector of maps from point to value, for each to_problem
467 : * Used for higher order elemental variables
468 : */
469 : void setSolutionVectorValues(const unsigned int var_index,
470 : const DofobjectToInterpValVec & dofobject_to_valsvec,
471 : const InterpCaches & interp_caches);
472 :
473 : /*
474 : * Modify the solution values after having set AND synchronized the solution
475 : * @param var the variable to set
476 : * @param dofobject_to_valsvec a vector of maps from DoF to values, for each to_problem
477 : * Used for nodal + constant monomial variables
478 : * @param interp_caches a vector of maps from point to value, for each to_problem
479 : * Used for higher order elemental variables
480 : */
481 : void correctSolutionVectorValues(const unsigned int var_index,
482 : const DofobjectToInterpValVec & dofobject_to_valsvec,
483 : const InterpCaches & interp_caches);
484 :
485 : /*
486 : * Cache pointInfo
487 : */
488 : void cacheOutgoingPointInfo(const Point point,
489 : const dof_id_type dof_object_id,
490 : const unsigned int problem_id,
491 : ProcessorToPointVec & outgoing_points);
492 :
493 : /**
494 : * Compute minimum distance
495 : * @param p the point of interest
496 : * @param bbox the bounding box to find the minimum distance from
497 : */
498 : Real bboxMinDistance(const Point & p, const BoundingBox & bbox) const;
499 :
500 : /**
501 : * Compute max distance
502 : * @param p the point of interest
503 : * @param bbox the bounding box to find the maximum distance from
504 : */
505 : Real bboxMaxDistance(const Point & p, const BoundingBox & bbox) const;
506 :
507 : /**
508 : * @brief Obtains the max dimensions to scale all points in the mesh
509 : * @return the maximum dimension in each coordinate axis of all target problems
510 : */
511 : Point getMaxToProblemsBBoxDimensions() const;
512 :
513 : /**
514 : * Get from bounding boxes for given domains and boundaries
515 : */
516 : std::vector<BoundingBox> getRestrictedFromBoundingBoxes() const;
517 :
518 : /**
519 : * Get global index for the first app each processes owns
520 : * Requires a global communication, must be called on every domain simultaneously
521 : */
522 : std::vector<unsigned int> getGlobalStartAppPerProc() const;
523 : };
524 :
525 : // Anonymous namespace for data, functors to use with GenericProjector.
526 : namespace GeneralFieldTransfer
527 : {
528 : // Transfer::OutOfMeshValue is an actual number. Why? Why!
529 : static_assert(std::numeric_limits<Real>::has_infinity,
530 : "What are you trying to use for Real? It lacks infinity!");
531 : extern Number OutOfMeshValue;
532 :
533 : inline bool
534 12360229 : isOutOfMeshValue(Number val)
535 : {
536 : // Might need to be changed for e.g. NaN
537 12360229 : return val == GeneralFieldTransfer::OutOfMeshValue;
538 : }
539 :
540 : // We need two functors that record point (value and gradient,
541 : // respectively) requests, so we know what queries we need to make
542 : // to other processors
543 :
544 : /**
545 : * Value request recording base class
546 : */
547 : template <typename Output>
548 : class RecordRequests
549 : {
550 : protected:
551 : typedef typename libMesh::TensorTools::MakeBaseNumber<Output>::type DofValueType;
552 :
553 : public:
554 : typedef typename libMesh::TensorTools::MakeReal<Output>::type RealType;
555 : typedef DofValueType ValuePushType;
556 : typedef Output FunctorValue;
557 :
558 808 : RecordRequests() {}
559 :
560 2091 : RecordRequests(RecordRequests & primary) : _primary(&primary) {}
561 :
562 2893 : ~RecordRequests()
563 : {
564 2893 : if (_primary)
565 : {
566 2091 : Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
567 4182 : _primary->_points_requested.insert(
568 2091 : _primary->_points_requested.end(), _points_requested.begin(), _points_requested.end());
569 2091 : }
570 2893 : }
571 :
572 2091 : void init_context(libMesh::FEMContext &) {}
573 :
574 0 : Output eval_at_node(const libMesh::FEMContext &,
575 : unsigned int /*variable_index*/,
576 : unsigned int /*elem_dim*/,
577 : const Node & n,
578 : bool /*extra_hanging_dofs*/,
579 : const Real /*time*/)
580 : {
581 0 : _points_requested.push_back(n);
582 0 : return 0;
583 : }
584 :
585 283200 : Output eval_at_point(const libMesh::FEMContext &,
586 : unsigned int /*variable_index*/,
587 : const Point & n,
588 : const Real /*time*/,
589 : bool /*skip_context_check*/)
590 : {
591 283200 : _points_requested.push_back(n);
592 283200 : return 0;
593 : }
594 :
595 53920 : bool is_grid_projection() { return false; }
596 :
597 0 : void eval_mixed_derivatives(const libMesh::FEMContext & /*c*/,
598 : unsigned int /*i*/,
599 : unsigned int /*dim*/,
600 : const Node & /*n*/,
601 : std::vector<Output> & /*derivs*/)
602 : {
603 0 : mooseError("Not implemented");
604 : } // this is only for grid projections
605 :
606 0 : void eval_old_dofs(
607 : const Elem &, unsigned int, unsigned int, std::vector<dof_id_type> &, std::vector<Output> &)
608 : {
609 0 : mooseError("Not implemented");
610 : }
611 :
612 0 : void eval_old_dofs(const Elem &,
613 : const libMesh::FEType &,
614 : unsigned int,
615 : unsigned int,
616 : std::vector<dof_id_type> &,
617 : std::vector<Output> &)
618 : {
619 0 : mooseError("Not implemented");
620 : }
621 :
622 802 : std::vector<Point> & points_requested() { return _points_requested; }
623 :
624 : private:
625 : /// Vector of points requested
626 : std::vector<Point> _points_requested;
627 :
628 : RecordRequests * _primary = nullptr;
629 : };
630 :
631 : // We need a null action functor to use
632 : // with them (because we won't be ready to set any values at that point)
633 : template <typename Val>
634 : class NullAction
635 : {
636 : public:
637 : typedef Val InsertInput;
638 :
639 404 : NullAction() {}
640 :
641 183456 : void insert(dof_id_type, Val) {}
642 :
643 0 : void insert(const std::vector<dof_id_type> &, const DenseVector<Val> &) {}
644 : };
645 :
646 : // We need two functors that respond to point (value and gradient,
647 : // respectively) requests based on the cached values of queries answered by
648 : // other processors.
649 :
650 : /**
651 : * Value request response base class
652 : */
653 : template <typename Output>
654 : class CachedData
655 : {
656 : protected:
657 : typedef typename libMesh::TensorTools::MakeBaseNumber<Output>::type DofValueType;
658 :
659 : public:
660 : typedef PointIndexedMap Cache;
661 :
662 : typedef typename libMesh::TensorTools::MakeReal<Output>::type RealType;
663 : typedef DofValueType ValuePushType;
664 : typedef Output FunctorValue;
665 :
666 : /**
667 : * Constructor
668 : * @param cache a map/cache to search for points in
669 : * @param backup a function that can be queried for a point value when the cache doesnt have it
670 : */
671 401 : CachedData(const Cache & cache, const libMesh::FunctionBase<Output> & backup, Real default_value)
672 401 : : _cache(cache), _backup(backup.clone()), _default_value(default_value)
673 : {
674 401 : }
675 :
676 : /// Copy constructor
677 2091 : CachedData(const CachedData & primary)
678 2091 : : _cache(primary._cache),
679 2091 : _backup(primary._backup->clone()),
680 2091 : _default_value(primary._default_value)
681 : {
682 2091 : }
683 :
684 2091 : void init_context(libMesh::FEMContext &) {}
685 :
686 : /// Gets a value at the node location
687 0 : Output eval_at_node(const libMesh::FEMContext &,
688 : unsigned int /*i*/,
689 : unsigned int /*elem_dim*/,
690 : const Node & n,
691 : bool /*extra_hanging_dofs*/,
692 : const Real /*time*/)
693 : {
694 0 : auto it = _cache.find(n);
695 0 : if (it == _cache.end())
696 : {
697 0 : if (_default_value != GeneralFieldTransfer::OutOfMeshValue)
698 0 : return _default_value;
699 : else
700 0 : return (*_backup)(n);
701 : }
702 : else
703 0 : return it->second;
704 : }
705 :
706 : /// Gets a value at a point
707 283200 : Output eval_at_point(const libMesh::FEMContext &,
708 : unsigned int /*i*/,
709 : const Point & n,
710 : const Real /*time*/,
711 : bool /*skip_context_check*/)
712 : {
713 283200 : auto it = _cache.find(n);
714 283200 : if (it == _cache.end())
715 : {
716 211152 : if (_default_value != GeneralFieldTransfer::OutOfMeshValue)
717 211152 : return _default_value;
718 : else
719 0 : return (*_backup)(n);
720 : }
721 : else
722 72048 : return it->second;
723 : }
724 :
725 53920 : bool is_grid_projection() { return false; }
726 :
727 0 : void eval_mixed_derivatives(const libMesh::FEMContext & /*c*/,
728 : unsigned int /*i*/,
729 : unsigned int /*dim*/,
730 : const Node & /*n*/,
731 : std::vector<Output> & /*derivs*/)
732 : {
733 0 : mooseError("Not implemented");
734 : } // this is only for grid projections
735 :
736 0 : void eval_old_dofs(
737 : const Elem &, unsigned int, unsigned int, std::vector<dof_id_type> &, std::vector<Output> &)
738 : {
739 0 : mooseError("Not implemented");
740 : }
741 :
742 0 : void eval_old_dofs(const Elem &,
743 : const libMesh::FEType &,
744 : unsigned int,
745 : unsigned int,
746 : std::vector<dof_id_type> &,
747 : std::vector<Output> &)
748 : {
749 0 : mooseError("Not implemented");
750 : }
751 :
752 : private:
753 : /// Data to return for cached points
754 : const Cache & _cache;
755 :
756 : /// Function to evaluate for uncached points
757 : std::unique_ptr<libMesh::FunctionBase<Output>> _backup;
758 :
759 : /// Default value when no point is found
760 : const Real _default_value;
761 : };
762 :
763 : }
|