https://mooseframework.inl.gov
Public Member Functions | Protected Attributes | List of all members
ThreadedRadialAverageLoop Class Reference

RadialAverage threaded loop. More...

#include <ThreadedRadialAverageLoop.h>

Public Member Functions

 ThreadedRadialAverageLoop (RadialAverage &)
 
 ThreadedRadialAverageLoop (const ThreadedRadialAverageLoop &x, Threads::split split)
 Splitting constructor. More...
 
virtual ~ThreadedRadialAverageLoop ()
 dummy virtual destructor More...
 
void operator() (const QPDataRange &range)
 parens operator with the code that is executed in threads More...
 
virtual void join (const ThreadedRadialAverageLoop &)
 thread join method More...
 

Protected Attributes

RadialAverage_radavg
 rasterizer to manage the sample data More...
 
THREAD_ID _tid
 ID number of the current thread. More...
 

Detailed Description

RadialAverage threaded loop.

Definition at line 22 of file ThreadedRadialAverageLoop.h.

Constructor & Destructor Documentation

◆ ThreadedRadialAverageLoop() [1/2]

ThreadedRadialAverageLoop::ThreadedRadialAverageLoop ( RadialAverage green)

Definition at line 25 of file ThreadedRadialAverageLoop.C.

25 : _radavg(green) {}
RadialAverage & _radavg
rasterizer to manage the sample data

◆ ThreadedRadialAverageLoop() [2/2]

ThreadedRadialAverageLoop::ThreadedRadialAverageLoop ( const ThreadedRadialAverageLoop x,
Threads::split  split 
)

Splitting constructor.

Definition at line 28 of file ThreadedRadialAverageLoop.C.

30  : _radavg(x._radavg)
31 {
32 }
RadialAverage & _radavg
rasterizer to manage the sample data

◆ ~ThreadedRadialAverageLoop()

virtual ThreadedRadialAverageLoop::~ThreadedRadialAverageLoop ( )
inlinevirtual

dummy virtual destructor

Definition at line 31 of file ThreadedRadialAverageLoop.h.

31 {}

Member Function Documentation

◆ join()

virtual void ThreadedRadialAverageLoop::join ( const ThreadedRadialAverageLoop )
inlinevirtual

thread join method

Definition at line 37 of file ThreadedRadialAverageLoop.h.

37 {}

◆ operator()()

void ThreadedRadialAverageLoop::operator() ( const QPDataRange range)

parens operator with the code that is executed in threads

Definition at line 35 of file ThreadedRadialAverageLoop.C.

36 {
37  // fetch data from parent
38  const auto radius = _radavg._radius;
39  const auto & qp_data = _radavg._qp_data;
40  const auto & kd_tree = _radavg._kd_tree;
41  const auto & weights_type = _radavg._weights_type;
42 
43  // tree search data structures
44  std::vector<nanoflann::ResultItem<std::size_t, Real>> ret_matches;
45  nanoflann::SearchParameters search_params;
46 
47  // result map entry
48  const auto end_it = _radavg._average.end();
49  auto it = end_it;
50 
51  // iterate over qp range
52  for (auto && local_qp : qpdata_range)
53  {
54  // Look up result map iterator only if we enter a new element. this saves a bunch
55  // of map lookups because same element entries are consecutive in the qp_data vector.
56  if (it == end_it || it->first != local_qp._elem_id)
57  it = _radavg._average.find(local_qp._elem_id);
58 
59  // initialize result entry
60  mooseAssert(it != end_it, "Current element id not found in result set.");
61  auto & sum = it->second[local_qp._qp];
62  sum = 0.0;
63 
64  ret_matches.clear();
65  std::size_t n_result =
66  kd_tree->radiusSearch(&(local_qp._q_point(0)), radius * radius, ret_matches, search_params);
67  Real total_vol = 0.0;
68  Real weight = 1.0;
69  for (std::size_t j = 0; j < n_result; ++j)
70  {
71  const auto & other_qp = qp_data[ret_matches[j].first];
72  switch (weights_type)
73  {
75  break;
76 
78  weight = radius - std::sqrt(ret_matches[j].second);
79  break;
80 
82  weight = std::cos(std::sqrt(ret_matches[j].second) / radius * libMesh::pi) + 1.0;
83  break;
84  }
85 
86  sum += other_qp._value * other_qp._volume * weight;
87  total_vol += other_qp._volume * weight;
88  }
89  sum /= total_vol;
90  }
91 }
const Real _radius
cut-off radius
Definition: RadialAverage.h:87
enum RadialAverage::WeightsType _weights_type
const Real radius
std::unique_ptr< KDTreeType > _kd_tree
spatial index (nanoflann guarantees this to be threadsafe under read-only operations) ...
RadialAverage & _radavg
rasterizer to manage the sample data
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template cos(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(cos
std::vector< QPData > _qp_data
gathered data
Definition: RadialAverage.h:93
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
Result _average
average result
Definition: RadialAverage.h:96
SearchParams SearchParameters
const Real pi

Member Data Documentation

◆ _radavg

RadialAverage& ThreadedRadialAverageLoop::_radavg
protected

rasterizer to manage the sample data

Definition at line 41 of file ThreadedRadialAverageLoop.h.

Referenced by operator()().

◆ _tid

THREAD_ID ThreadedRadialAverageLoop::_tid
protected

ID number of the current thread.

Definition at line 44 of file ThreadedRadialAverageLoop.h.


The documentation for this class was generated from the following files: