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);
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 991 : 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 :
174 : /**
175 : * Whether a point lies inside the mesh division delineated by the MeshDivision object
176 : * @param pt point to examine, in the local coordinates (source frame for from_direction=true)
177 : * @param i_local the index of the problem to consider, holding the mesh division to examine
178 : * @param only_from_this_mesh_div a mesh division index that must be matched when the
179 : * to/from_mesh_division_behavior for the direction examined is MATCH_DIVISION/SUBAPP_INDEX
180 : * It is ignored otherwise
181 : */
182 : bool acceptPointMeshDivision(const Point & pt,
183 : const unsigned int i_local,
184 : const unsigned int only_from_this_mesh_div) const;
185 :
186 : /**
187 : * Whether a point is closest to a position at the index specified than any other position
188 : * @param pos_index the index of the position to consider in the positions vector
189 : * @param pt the point
190 : * @return whether the point is closest to this position than any other in the positions vector
191 : */
192 : bool closestToPosition(unsigned int pos_index, const Point & pt) const;
193 :
194 : /// Origin array/vector variable components
195 : const std::vector<unsigned int> _from_var_components;
196 :
197 : /// Target array/vector variable components
198 : const std::vector<unsigned int> _to_var_components;
199 :
200 : /// Whether to use bounding boxes to determine the applications that may receive point requests
201 : /// then send value data, and at other various checks
202 : const bool _use_bounding_boxes;
203 :
204 : /// Whether to keep track of the distance from the requested point to the app position
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
215 : const Positions * _nearest_positions_obj;
216 :
217 : /// Whether the source app mesh must actually contain the points for them to be considered or whether
218 : /// the bounding box is enough. If false, we can interpolate between apps
219 : bool _source_app_must_contain_point;
220 :
221 : /// Origin block(s) restriction
222 : std::set<SubdomainID> _from_blocks;
223 :
224 : /// Target block(s) restriction
225 : std::set<SubdomainID> _to_blocks;
226 :
227 : /// Target boundary(ies) restriction
228 : std::set<BoundaryID> _to_boundaries;
229 :
230 : /// Origin boundary(ies) restriction
231 : std::set<BoundaryID> _from_boundaries;
232 :
233 : /// Division of the origin mesh
234 : std::vector<const MeshDivision *> _from_mesh_divisions;
235 :
236 : /// Division of the target mesh
237 : std::vector<const MeshDivision *> _to_mesh_divisions;
238 :
239 : /// How to use the origin mesh divisions to restrict the transfer
240 : const MooseEnum & _from_mesh_division_behavior;
241 :
242 : /// How to use the target mesh divisions to restrict the transfer
243 : const MooseEnum & _to_mesh_division_behavior;
244 :
245 : /// Matching enum for the mesh division behaviors.
246 : enum MeshDivisionTransferUse
247 : {
248 : RESTRICTION,
249 : MATCH_DIVISION_INDEX,
250 : MATCH_SUBAPP_INDEX
251 : };
252 :
253 : /// Whether elemental variable boundary restriction is considered by element side or element nodes
254 : const bool _elemental_boundary_restriction_on_sides;
255 :
256 : /// Point locators, useful to examine point location with regards to domain restriction
257 : std::vector<std::unique_ptr<libMesh::PointLocatorBase>> _from_point_locators;
258 :
259 : /// First app each processor owns, indexed by processor
260 : /// If no app on the processor, will have a -1 for the app start instead
261 : std::vector<unsigned int> _global_app_start_per_proc;
262 :
263 : /// Whether or not a greedy strategy will be used
264 : /// If true, all the partitions will be checked for a given
265 : /// outgoing point
266 : bool _greedy_search;
267 :
268 : /// Whether to look for conflicts between origin points, multiple valid values for a target point
269 : bool _search_value_conflicts;
270 :
271 : /// Whether we already output the search value conflicts
272 : bool _already_output_search_value_conflicts;
273 :
274 : /// How many conflicts are output to console
275 : const unsigned int _search_value_conflicts_max_log;
276 :
277 : /**
278 : * @brief Detects whether two source values are valid and equidistant for a desired target
279 : * location
280 : * @param value_1 value from the first value source / subapp
281 : * @param value_2 value from the second value source / subapp
282 : * @param distance_1 distance from the first source
283 : * @param distance_2 distance from the second source
284 : * @return true if the values are different and distances from the source points/apps are the same
285 : */
286 : bool detectConflict(Real value_1, Real value_2, Real distance_1, Real distance_2) const;
287 :
288 : /**
289 : * Register a potential value conflict, e.g. two or more equidistant source points for a single
290 : * target point, with different values possible
291 : * @param problem problem ID for the point of interest.
292 : * For local conflicts, use origin problem id, for received conflicts, use target id
293 : * @param dof_id id id of the DoF is transferring a DoF. If not, use -1
294 : * @param p point where the conflict happens
295 : * @param dist distance between the origin and the target
296 : * @param local if true, local conflict found when gathering data to send, if false,
297 : * received value conflict found when receiving data from multiple source problems
298 : */
299 : void registerConflict(unsigned int problem, dof_id_type dof_id, Point p, Real dist, bool local);
300 :
301 : private:
302 : /// The target variables
303 : std::vector<MooseVariableFieldBase *> _to_variables;
304 :
305 : /// A map from pid to a set of points
306 : typedef std::unordered_map<processor_id_type, std::vector<std::pair<Point, unsigned int>>>
307 : ProcessorToPointVec;
308 :
309 : /// Point information
310 : struct PointInfo
311 : {
312 : unsigned int problem_id; // problem id
313 : dof_id_type dof_object_id; // node or elem id
314 : };
315 :
316 : /// InterpInfo
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 :
324 : /// A map from pid to a set of point info
325 : typedef std::unordered_map<processor_id_type, std::vector<PointInfo>> ProcessorToPointInfoVec;
326 :
327 : /// A vector, indexed by to-problem id, of maps from dof object to interpolation values
328 : typedef std::vector<std::unordered_map<dof_id_type, InterpInfo>> DofobjectToInterpValVec;
329 :
330 : /// A map from Point to interpolation values
331 : /// NOTE: this is not an asynchronous cache. It is built to completion during the transfer
332 : /// and used as a whole to reconstruct the target variable
333 : typedef PointIndexedMap InterpCache;
334 :
335 : /// A vector of such caches, indexed by to_problem
336 : typedef std::vector<InterpCache> InterpCaches;
337 :
338 : /// The number of variables to transfer
339 : unsigned int _var_size;
340 :
341 : /// Error out when some points can not be located
342 : bool _error_on_miss;
343 :
344 : /// Value to use when no received data is valid for a target location
345 : const Real _default_extrapolation_value;
346 :
347 : /// How much we should relax bounding boxes
348 : Real _bbox_factor;
349 :
350 : /// Set the bounding box sizes manually
351 : std::vector<Real> _fixed_bbox_size;
352 :
353 : /// Number of source/from applications per processor. This vector is indexed by processor id
354 : std::vector<unsigned int> _froms_per_proc;
355 :
356 : /// Bounding boxes for all source applications. The indexing of this vector is similar to
357 : ///'processor_id * n_local_subapps + i_local_subapp', except the number of local subapps can
358 : /// be different on each processor. Use the _from_per_proc vector to find the start/end index for a given
359 : /// processor. See MultiAppGeneralFieldTransfer::locatePointReceivers() for an example
360 : std::vector<BoundingBox> _from_bboxes;
361 :
362 : /// A map from processor to pointInfo vector
363 : ProcessorToPointInfoVec _processor_to_pointInfoVec;
364 :
365 : /// Keeps track of all local equidistant points to requested points, creating an indetermination
366 : /// in which values should be sent for that request
367 : /// We keep the origin problem ID, the dof ID, the point, and the distance origin-target
368 : /// If using nearest-positions the origin problem ID is not set
369 : std::vector<std::tuple<unsigned int, dof_id_type, Point, Real>> _local_conflicts;
370 :
371 : /// Keeps track of all received conflicts. Multiple problems (different subapps for example)
372 : /// are sending values for a target point that do not match and are equally valid/distant
373 : /// We keep the target problem ID, the point/dof ID, the point, and the origin-target distance.
374 : /// The distance indicates whether a potential conflict ended up materializing
375 : std::vector<std::tuple<unsigned int, dof_id_type, Point, Real>> _received_conflicts;
376 :
377 : /**
378 : * Initialize supporting attributes like bounding boxes, processor app indexes etc
379 : */
380 : void prepareToTransfer();
381 :
382 : /**
383 : * Performs the transfer for the variable of index i
384 : */
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 : */
400 : void cacheIncomingInterpVals(
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 :
410 : /**
411 : * Remove potential value conflicts that did not materialize because another source was closer
412 : * Several equidistant valid values were received, but they were not closest
413 : * @param var_index the index of the variable of interest
414 : * @param dofobject_to_valsvec a data structure mapping dofobjects to received values and
415 : * distances (used for nodal-value-dof-only variables and constant monomials)
416 : * @param distance_caches a cache holding the distances received (used for higher order elemental
417 : * variables)
418 : */
419 : void examineReceivedValueConflicts(const unsigned int var_index,
420 : const DofobjectToInterpValVec & dofobject_to_valsvec,
421 : const InterpCaches & distance_caches);
422 :
423 : /**
424 : * Remove potential value conflicts that did not materialize because another source was closer
425 : * Several equidistant valid values were found when computing values to send, but they were not
426 : * closest, another value got selected
427 : * @param var_index the index of the variable of interest
428 : * @param dofobject_to_valsvec a data structure mapping dofobjects to received values and
429 : * distances (used for nodal-value-dof-only variables and constant monomials)
430 : * @param distance_caches a cache holding the distances received (used for higher order elemental
431 : * variables)
432 : */
433 : void examineLocalValueConflicts(const unsigned int var_index,
434 : const DofobjectToInterpValVec & dofobject_to_valsvec,
435 : const InterpCaches & distance_caches);
436 :
437 : /// Report on conflicts between overlapping child apps, equidistant origin points etc
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 :
462 : /**
463 : * Compute minimum distance
464 : * @param p the point of interest
465 : * @param bbox the bounding box to find the minimum distance from
466 : */
467 : Real bboxMinDistance(const Point & p, const BoundingBox & bbox) const;
468 :
469 : /**
470 : * Compute max distance
471 : * @param p the point of interest
472 : * @param bbox the bounding box to find the maximum distance from
473 : */
474 : Real bboxMaxDistance(const Point & p, const BoundingBox & bbox) const;
475 :
476 : /**
477 : * @brief Obtains the max dimensions to scale all points in the mesh
478 : * @return the maximum dimension in each coordinate axis of all target problems
479 : */
480 : Point getMaxToProblemsBBoxDimensions() const;
481 :
482 : /**
483 : * Get from bounding boxes for given domains and boundaries
484 : */
485 : std::vector<BoundingBox> getRestrictedFromBoundingBoxes() const;
486 :
487 : /**
488 : * Get global index for the first app each processes owns
489 : * Requires a global communication, must be called on every domain simultaneously
490 : */
491 : std::vector<unsigned int> getGlobalStartAppPerProc() const;
492 : };
493 :
494 : // Anonymous namespace for data, functors to use with GenericProjector.
495 : namespace GeneralFieldTransfer
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!");
500 : extern Number BetterOutOfMeshValue;
501 :
502 : inline bool
503 12629882 : isBetterOutOfMeshValue(Number val)
504 : {
505 : // Might need to be changed for e.g. NaN
506 12629882 : return val == GeneralFieldTransfer::BetterOutOfMeshValue;
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 :
513 : /**
514 : * Value request recording base class
515 : */
516 : template <typename Output>
517 : class RecordRequests
518 : {
519 : protected:
520 : typedef typename libMesh::TensorTools::MakeBaseNumber<Output>::type DofValueType;
521 :
522 : public:
523 : typedef typename libMesh::TensorTools::MakeReal<Output>::type RealType;
524 : typedef DofValueType ValuePushType;
525 : typedef Output FunctorValue;
526 :
527 670 : RecordRequests() {}
528 :
529 1725 : RecordRequests(RecordRequests & primary) : _primary(&primary) {}
530 :
531 2387 : ~RecordRequests()
532 : {
533 2387 : if (_primary)
534 : {
535 1725 : Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
536 3450 : _primary->_points_requested.insert(
537 1725 : _primary->_points_requested.end(), _points_requested.begin(), _points_requested.end());
538 1725 : }
539 2387 : }
540 :
541 1725 : void init_context(libMesh::FEMContext &) {}
542 :
543 0 : Output eval_at_node(const libMesh::FEMContext &,
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 0 : _points_requested.push_back(n);
551 0 : return 0;
552 : }
553 :
554 265200 : Output eval_at_point(const libMesh::FEMContext &,
555 : unsigned int /*variable_index*/,
556 : const Point & n,
557 : const Real /*time*/,
558 : bool /*skip_context_check*/)
559 : {
560 265200 : _points_requested.push_back(n);
561 265200 : return 0;
562 : }
563 :
564 47040 : bool is_grid_projection() { return false; }
565 :
566 0 : void eval_mixed_derivatives(const libMesh::FEMContext & /*c*/,
567 : unsigned int /*i*/,
568 : unsigned int /*dim*/,
569 : const Node & /*n*/,
570 : std::vector<Output> & /*derivs*/)
571 : {
572 0 : mooseError("Not implemented");
573 : } // this is only for grid projections
574 :
575 0 : void eval_old_dofs(
576 : const Elem &, unsigned int, unsigned int, std::vector<dof_id_type> &, std::vector<Output> &)
577 : {
578 0 : mooseError("Not implemented");
579 : }
580 :
581 0 : 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 0 : mooseError("Not implemented");
589 : }
590 :
591 662 : std::vector<Point> & points_requested() { return _points_requested; }
592 :
593 : private:
594 : /// Vector of points requested
595 : std::vector<Point> _points_requested;
596 :
597 : RecordRequests * _primary = nullptr;
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>
603 : class NullAction
604 : {
605 : public:
606 : typedef Val InsertInput;
607 :
608 335 : NullAction() {}
609 :
610 170592 : void insert(dof_id_type, Val) {}
611 :
612 0 : 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 :
619 : /**
620 : * Value request response base class
621 : */
622 : template <typename Output>
623 : class CachedData
624 : {
625 : protected:
626 : typedef typename libMesh::TensorTools::MakeBaseNumber<Output>::type DofValueType;
627 :
628 : public:
629 : typedef PointIndexedMap Cache;
630 :
631 : typedef typename libMesh::TensorTools::MakeReal<Output>::type RealType;
632 : typedef DofValueType ValuePushType;
633 : typedef Output FunctorValue;
634 :
635 : /**
636 : * Constructor
637 : * @param cache a map/cache to search for points in
638 : * @param backup a function that can be queried for a point value when the cache doesnt have it
639 : */
640 331 : CachedData(const Cache & cache, const libMesh::FunctionBase<Output> & backup, Real default_value)
641 331 : : _cache(cache), _backup(backup.clone()), _default_value(default_value)
642 : {
643 331 : }
644 :
645 : /// Copy constructor
646 1725 : CachedData(const CachedData & primary)
647 1725 : : _cache(primary._cache),
648 1725 : _backup(primary._backup->clone()),
649 1725 : _default_value(primary._default_value)
650 : {
651 1725 : }
652 :
653 1725 : void init_context(libMesh::FEMContext &) {}
654 :
655 : /// Gets a value at the node location
656 0 : Output eval_at_node(const libMesh::FEMContext &,
657 : unsigned int /*i*/,
658 : unsigned int /*elem_dim*/,
659 : const Node & n,
660 : bool /*extra_hanging_dofs*/,
661 : const Real /*time*/)
662 : {
663 0 : auto it = _cache.find(n);
664 0 : if (it == _cache.end())
665 : {
666 0 : if (_default_value != GeneralFieldTransfer::BetterOutOfMeshValue)
667 0 : return _default_value;
668 : else
669 0 : return (*_backup)(n);
670 : }
671 : else
672 0 : return it->second;
673 : }
674 :
675 : /// Gets a value at a point
676 265200 : Output eval_at_point(const libMesh::FEMContext &,
677 : unsigned int /*i*/,
678 : const Point & n,
679 : const Real /*time*/,
680 : bool /*skip_context_check*/)
681 : {
682 265200 : auto it = _cache.find(n);
683 265200 : if (it == _cache.end())
684 : {
685 204128 : if (_default_value != GeneralFieldTransfer::BetterOutOfMeshValue)
686 204128 : return _default_value;
687 : else
688 0 : return (*_backup)(n);
689 : }
690 : else
691 61072 : return it->second;
692 : }
693 :
694 47040 : bool is_grid_projection() { return false; }
695 :
696 0 : void eval_mixed_derivatives(const libMesh::FEMContext & /*c*/,
697 : unsigned int /*i*/,
698 : unsigned int /*dim*/,
699 : const Node & /*n*/,
700 : std::vector<Output> & /*derivs*/)
701 : {
702 0 : mooseError("Not implemented");
703 : } // this is only for grid projections
704 :
705 0 : void eval_old_dofs(
706 : const Elem &, unsigned int, unsigned int, std::vector<dof_id_type> &, std::vector<Output> &)
707 : {
708 0 : mooseError("Not implemented");
709 : }
710 :
711 0 : 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 0 : mooseError("Not implemented");
719 : }
720 :
721 : private:
722 : /// Data to return for cached points
723 : const Cache & _cache;
724 :
725 : /// Function to evaluate for uncached points
726 : std::unique_ptr<libMesh::FunctionBase<Output>> _backup;
727 :
728 : /// Default value when no point is found
729 : const Real _default_value;
730 : };
731 :
732 : }
|