https://mooseframework.inl.gov
MultiAppGeneralFieldNearestLocationTransfer.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"
15 
21 
22 {
23 public:
25 
27 
28  void initialSetup() override;
29 
30  // Use solution invalid output for these warnings
32 
33 protected:
34  virtual void prepareEvaluationOfInterpValues(const unsigned int var_index) override;
35 
36  virtual void
37  evaluateInterpValues(const std::vector<std::pair<Point, unsigned int>> & incoming_points,
38  std::vector<std::pair<Real, Real>> & outgoing_vals) override;
39 
41  bool inBlocks(const std::set<SubdomainID> & blocks,
42  const MooseMesh & mesh,
43  const Elem * elem) const override;
44 
45 private:
46  bool usesMooseAppCoordTransform() const override { return true; }
47  /*
48  * Build KD-Trees for each local app
49  * @param var_index the index of the variable being transferred
50  * @details fills _local_kdtrees, _local_points and _local_values
51  * Indexing is: local apps (outer-indexing) OR positions (if using nearest_positions),
52  * local nodes (inner-indexing)
53  */
54  void buildKDTrees(const unsigned int var_index);
55 
56  /*
57  * Evaluate interpolation values for incoming points
58  * @param incoming_points all the points at which we need values
59  * @param outgoing_vals vector containing the values and distances from point to nearest node
60  */
62  const std::vector<std::pair<Point, unsigned int>> & incoming_points,
63  std::vector<std::pair<Real, Real>> & outgoing_vals);
64 
67  void computeNumSources();
68 
70  unsigned int getAppIndex(unsigned int kdtree_index, unsigned int app_index) const;
71 
73  unsigned int getNumAppsPerTree() const;
74 
76  unsigned int getNumDivisions() const;
77 
79  Point getPointInLocalSourceFrame(unsigned int i_from, const Point & pt) const;
80 
89  bool checkRestrictionsForSource(const Point & pt,
90  const unsigned int valid_mesh_div,
91  const unsigned int i_from) const;
92 
94  std::vector<std::shared_ptr<KDTree>> _local_kdtrees;
95 
97  std::vector<std::vector<std::shared_ptr<KDTree>>> _local_positions_kdtrees;
98 
100  std::vector<std::vector<Point>> _local_points;
101 
103  std::vector<std::vector<Real>> _local_values;
104 
106  unsigned int _num_nearest_points;
107 
109  const bool _group_subapps;
110 
112  std::vector<bool> _source_is_nodes;
113 
115  std::vector<bool> _use_zero_dof_for_value;
116 
118  unsigned int _num_sources;
119 };
bool inBlocks(const std::set< SubdomainID > &blocks, const MooseMesh &mesh, const Elem *elem) const override
const bool _group_subapps
Whether to group data when creating the nearest-point regions.
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.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
Point getPointInLocalSourceFrame(unsigned int i_from, const Point &pt) const
Transform a point towards the local frame.
std::vector< bool > _source_is_nodes
Whether the source of the values is at nodes (true) or centroids (false) for each variable...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
unsigned int getNumAppsPerTree() const
Number of applications which contributed nearest-locations to each KD-tree.
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:92
bool usesMooseAppCoordTransform() const override
Whether this transfer handles non-translation-based transformations, e.g.
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.
std::vector< std::shared_ptr< KDTree > > _local_kdtrees
KD-Trees for all the local source apps.
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...
virtual void evaluateInterpValues(const std::vector< std::pair< Point, unsigned int >> &incoming_points, std::vector< std::pair< Real, Real >> &outgoing_vals) override
std::vector< std::vector< std::shared_ptr< KDTree > > > _local_positions_kdtrees
KD-Trees for nodes nearest to a given position on each local source app.
virtual void prepareEvaluationOfInterpValues(const unsigned int var_index) override
Performs a geometric interpolation based on the values at the nearest nodes to a target location in t...
unsigned int getNumDivisions() const
Number of divisions (nearest-positions or source mesh divisions) used when building KD-Trees...
void initialSetup() override
Method called at the beginning of the simulation for checking integrity or doing one-time setup...
void computeNumSources()
Pre-compute the number of sources Number of KDTrees used to hold the locations and variable value dat...
It is a general field transfer.
bool inBlocks(const std::set< SubdomainID > &blocks, const Elem *elem) const