https://mooseframework.inl.gov
RadialAverage.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 
12 #include "ElementUserObject.h"
13 #include "DataIO.h"
14 #include "PointListAdaptor.h"
15 #include "Function.h"
16 
17 #include "libmesh/data_type.h"
18 
19 #include <set>
20 #include <map>
21 #include <vector>
22 #include <memory>
23 
25 
34 {
35 public:
37 
39 
40  virtual void initialSetup() override;
41 
42  virtual void initialize() override;
43  virtual void execute() override;
44  virtual void finalize() override;
45  virtual void threadJoin(const UserObject & y) override;
46  virtual void meshChanged() override;
47 
49  struct QPData
50  {
52  Point _q_point;
56  short _qp;
61 
63  QPData(const Point & q_point, dof_id_type elem_id, short qp, Real volume, Real value)
64  : _q_point(q_point), _elem_id(elem_id), _qp(qp), _volume(volume), _value(value)
65  {
66  }
67  };
68 
69  using Result = std::map<dof_id_type, std::vector<Real>>;
70  const Result & getAverage() const { return _average; }
71 
72 protected:
74 
76  enum class WeightsType
77  {
78  CONSTANT,
79  LINEAR,
80  COSINE
81  } _weights_type;
82 
85 
87  const Real _radius;
88 
90  const Real _padding;
91 
93  std::vector<QPData> _qp_data;
94 
97 
98  using KDTreeType = nanoflann::KDTreeSingleIndexAdaptor<
99  nanoflann::L2_Simple_Adaptor<Real, PointListAdaptor<QPData>>,
101  LIBMESH_DIM,
102  std::size_t>;
103 
105  std::unique_ptr<KDTreeType> _kd_tree;
106 
108  std::vector<std::vector<const Elem *>> _nodes_to_elem_map;
109 
111  std::set<Point> _boundary_nodes;
112 
114  std::set<std::size_t> _boundary_data_indices;
115 
117  std::vector<std::set<std::size_t>> _communication_lists;
119 
121  std::vector<processor_id_type> _candidate_procs;
122 
124 
130 
132 };
133 
134 namespace TIMPI
135 {
136 
137 template <>
138 class StandardType<RadialAverage::QPData> : public DataType
139 {
140 public:
141  explicit StandardType(const RadialAverage::QPData * example = nullptr);
143  ~StandardType() { this->free(); }
144 };
145 
146 } // namespace TIMPI
const Real _radius
cut-off radius
Definition: RadialAverage.h:87
virtual void initialSetup() override
Gets called at the beginning of the simulation before this object is asked to do its job...
Definition: RadialAverage.C:87
Real _volume
current value * _JxW
Definition: RadialAverage.h:58
enum RadialAverage::WeightsType _weights_type
std::set< Point > _boundary_nodes
set of nodes on the boundary of the current processor domain
std::vector< std::set< std::size_t > > _communication_lists
QPData indices to send to the various processors.
const MaterialProperty< Real > & _prop
material to be averaged
Definition: RadialAverage.h:84
const unsigned int invalid_uint
const Result & getAverage() const
Definition: RadialAverage.h:70
std::unique_ptr< KDTreeType > _kd_tree
spatial index (nanoflann guarantees this to be threadsafe under read-only operations) ...
std::map< dof_id_type, std::vector< Real > > Result
Definition: RadialAverage.h:69
virtual void finalize() override
Finalize.
virtual void execute() override
Execute method.
virtual void threadJoin(const UserObject &y) override
Must override.
quadrature point data
Definition: RadialAverage.h:49
WeightsType
distance based weight function
Definition: RadialAverage.h:76
QPData(const Point &q_point, dof_id_type elem_id, short qp, Real volume, Real value)
Definition: RadialAverage.h:63
Real _value
variable value
Definition: RadialAverage.h:60
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
std::set< std::size_t > _boundary_data_indices
set of all _qp_data indices that are within _radius of any _boundary_nodes
processor_id_type _my_pid
void updateCommunicationLists()
virtual void initialize() override
Called before execute() is ever called so that data can be cleared.
Definition: RadialAverage.C:94
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
const Real _padding
communication padding distance
Definition: RadialAverage.h:90
StandardType(const T *example=nullptr)
std::vector< std::vector< const Elem * > > _nodes_to_elem_map
The data structure used to find neighboring elements give a node ID.
static InputParameters validParams()
Definition: RadialAverage.C:45
short _qp
index of the quadrature point
Definition: RadialAverage.h:56
unsigned int PerfID
Definition: MooseTypes.h:212
nanoflann::KDTreeSingleIndexAdaptor< nanoflann::L2_Simple_Adaptor< Real, PointListAdaptor< QPData > >, PointListAdaptor< QPData >, LIBMESH_DIM, std::size_t > KDTreeType
virtual void meshChanged() override
Called on this object when the mesh changes.
uint8_t processor_id_type
RadialAverage(const InputParameters &parameters)
Definition: RadialAverage.C:72
PerfID _perf_updatelists
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
std::vector< QPData > _qp_data
gathered data
Definition: RadialAverage.h:93
bool _update_communication_lists
PerfID _perf_meshchanged
RadialAverage threaded loop.
Real volume(const MeshBase &mesh, unsigned int dim=libMesh::invalid_uint)
Point _q_point
physical coordinates of the quadrature point
Definition: RadialAverage.h:52
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Result _average
average result
Definition: RadialAverage.h:96
dof_id_type _elem_id
element id
Definition: RadialAverage.h:54
Gather and communicate a full list of all quadrature points and the values of a selected material pro...
Definition: RadialAverage.h:33
const InputParameters & parameters() const
Get the parameters of the object.
PerfID _perf_finalize
Base class for user-specific data.
Definition: UserObject.h:40
std::vector< processor_id_type > _candidate_procs
processors to send (potentially empty) data to
uint8_t dof_id_type