OGS
CellAverageAlgorithm.h
Go to the documentation of this file.
1
11#pragma once
12
13#include "CellAverageData.h"
15
16namespace ProcessLib
17{
18namespace detail
19{
20void computeCellAverages(CellAverageData& cell_average_data,
21 std::string const& name, unsigned const num_comp,
22 auto&& flattened_ip_data_accessor,
23 auto const& local_assemblers)
24{
25 auto& prop_vec =
26 cell_average_data.getOrCreatePropertyVector(name, num_comp);
27
28 for (std::size_t i = 0; i < local_assemblers.size(); ++i)
29 {
30 auto const& loc_asm = *local_assemblers[i];
31 auto const& ip_data = flattened_ip_data_accessor(loc_asm);
32 assert(ip_data.size() % num_comp == 0);
33 auto const num_ips =
34 static_cast<Eigen::Index>(ip_data.size() / num_comp);
35 Eigen::Map<const Eigen::MatrixXd> ip_data_mapped{ip_data.data(),
36 num_comp, num_ips};
37
38 Eigen::Map<Eigen::VectorXd>{&prop_vec[i * num_comp], num_comp} =
39 ip_data_mapped.rowwise().mean();
40 }
41}
42} // namespace detail
43
44template <int dim, typename LAIntf>
46 CellAverageData& cell_average_data,
47 std::vector<std::unique_ptr<LAIntf>> const& local_assemblers)
48{
49 auto const callback = [&cell_average_data, &local_assemblers](
50 std::string const& name,
51 unsigned const num_comp,
52 auto&& flattened_ip_data_accessor)
53 {
54 detail::computeCellAverages(cell_average_data, name, num_comp,
55 flattened_ip_data_accessor,
56 local_assemblers);
57 };
58
60 LAIntf>(
61 LAIntf::getReflectionDataForOutput(), callback);
62}
63} // namespace ProcessLib
void forEachReflectedFlattenedIPDataAccessor(ReflData const &reflection_data, Callback const &callback)
void computeCellAverages(CellAverageData &cell_average_data, std::string const &name, unsigned const num_comp, auto &&flattened_ip_data_accessor, auto const &local_assemblers)
void computeCellAverages(CellAverageData &cell_average_data, std::vector< std::unique_ptr< LAIntf > > const &local_assemblers)
MeshLib::PropertyVector< double > & getOrCreatePropertyVector(std::string const &name, unsigned const num_comp)