MaMiCo 1.2
Loading...
Searching...
No Matches
FilterInterface.h
1// This file is part of the Mamico project. For conditions of distribution
2// and use, please see the copyright notice in Mamico's main folder, or at
3// www5.in.tum.de/mamico
4
5#pragma once
6#include <array>
7#include <vector>
8
9#include "coupling/indexing/CellIndex.h"
10
11// #define DEBUG_FILTER_INTERFACE
12
13namespace coupling {
14namespace filtering {
15template <unsigned int dim> class FilterInterface;
16}
17} // namespace coupling
18
32template <unsigned int dim> class coupling::filtering::FilterInterface {
33public:
34 /*
35 * Filter constructors are called during instanciation of their corresponding
36 * FilterSequence. You can customize parameterization in
37 * coupling::FilterSequence::loadFiltersFromXML(...).
38 */
39 FilterInterface(const std::vector<coupling::datastructures::CouplingCell<dim>*>& inputCellVector,
40 const std::vector<coupling::datastructures::CouplingCell<dim>*>& outputCellVector, const std::array<bool, 7> filteredValues, const char* type)
41 :
42
43 _inputCells(inputCellVector), _outputCells(outputCellVector), _type(type) {
44 // microscopic mass
45 if (filteredValues[0]) {
46 _scalarAccessFunctionPairs.push_back(
48 }
49 // microscopic momentum
50 if (filteredValues[1]) {
51 _vectorAccessFunctionPairs.push_back(
53 }
54 // macroscopic mass
55 if (filteredValues[2]) {
56 _scalarAccessFunctionPairs.push_back(
58 }
59 // macroscopic momentum
60 if (filteredValues[3]) {
61 _vectorAccessFunctionPairs.push_back(
63 }
64 // potential energy
65 if (filteredValues[4]) {
66 _scalarAccessFunctionPairs.push_back(
68 }
69 // velocity
70 if (filteredValues[5]) {
71 _vectorAccessFunctionPairs.push_back(
73 }
74 // temperature
75 if (filteredValues[6]) {
76 _scalarAccessFunctionPairs.push_back(
78 }
79 }
80
81 FilterInterface(const char* type)
82 : _type(type) { /* Used by incomplete implementations of FilterInterface.
83 Should be redesigned via meta class.*/
84 }
85
86 virtual ~FilterInterface() {};
87
88 /*
89 * Applies the filter to all cells that are within the filter's sequence's
90 * domain.
91 *
92 * It is very important that this method provides complete output data,
93 * i.e uses all elements of _scalarSetters and _vectorSetters on all elements
94 * of _outputCells. If this is not the case, you dont want to use this
95 * interface, but rather coupling::FilterInterfaceReadOnly and use its method
96 * copyInputToOutput().
97 */
98 virtual void operator()() = 0;
99
100 void updateCellData(const std::vector<coupling::datastructures::CouplingCell<dim>*>& new_inputCells,
101 const std::vector<coupling::datastructures::CouplingCell<dim>*>& new_outputCells) {
102 if (new_inputCells.size() != new_outputCells.size())
103 throw std::runtime_error("New input-, output-, and indexing vectors must "
104 "be of identical size.");
105
106 _inputCells = new_inputCells;
107 _outputCells = new_outputCells;
108
109#ifdef DEBUG_FILTER_INTERFACE
110 std::cout << " FI: Updated cell data." << std::endl;
111#endif
112 }
113
114 /*
115 * Basic Getters/Setters
116 */
117 const char* getType() const { return _type; }
118 std::vector<coupling::datastructures::CouplingCell<dim>*> getInputCells() const { return _inputCells; }
119 std::vector<coupling::datastructures::CouplingCell<dim>*> getOutputCells() const { return _outputCells; }
120
121 using CellIndex_T = coupling::indexing::CellIndex<dim, coupling::indexing::IndexTrait::local, coupling::indexing::IndexTrait::md2macro,
122 coupling::indexing::IndexTrait::noGhost>;
123 /*
124 * Advanced Getters/Setters
125 */
126 coupling::datastructures::CouplingCell<dim>* getInputCellOfIndex(const CellIndex_T& index) {
127 if (index.get() < _inputCells.size) {
128 return _inputCells[index.get()];
129 } else {
130 std::cout << "Index not found: " << index << std::endl;
131 throw std::runtime_error("FilterInterface: getInputCellofIndex(): Could not find index.");
132 }
133 }
134 coupling::datastructures::CouplingCell<dim>* getOutputCellOfIndex(const CellIndex_T& index) {
135 if (index.get() < _outputCells.size) {
136 return _outputCells[index.get()];
137 } else {
138 std::cout << "Index not found: " << index << std::endl;
139 throw std::runtime_error("FilterInterface: getOutputCellofIndex(): Could not find index.");
140 }
141 }
142
143 /*
144 * Only used in one scenario:
145 * - this is at index 0 in FS
146 * - new filter gets dynamically linked into FS at index 0
147 * In that case, this was previously getting input from MD but won't be any
148 * longer. The newly added filter will provide input for this one instead.
149 */
150 void setInputCells(const std::vector<coupling::datastructures::CouplingCell<dim>*>& newInputCells) { _inputCells = newInputCells; }
151
152 // Size = number of cells in this filter.
153 int getSize() const { return _inputCells.size(); }
154
155 /*
156 * Used by filter implementations to iterate over physical properties stored
157 * in a CouplingCell.
158 *
159 * Examplary usage:
160 * 'for(auto scalar : _scalarAccessFunctionPairs)' loops over all scalar
161 * properties (e.g. macro/micro mass, temperate) filtered by this specific
162 * filter.
163 */
165 const double& (coupling::datastructures::CouplingCell<dim>::*get)() const; // getter function pointer
166 void (coupling::datastructures::CouplingCell<dim>::*set)(const double&); // setter function pointer
167 };
172
173protected:
181 std::vector<coupling::datastructures::CouplingCell<dim>*> _inputCells;
182 std::vector<coupling::datastructures::CouplingCell<dim>*> _outputCells;
183
184 // scalars getters/setters
185 std::vector<ScalarAccessFunctionPair> _scalarAccessFunctionPairs;
186
187 // vectors getters/setters
188 std::vector<VectorAccessFunctionPair> _vectorAccessFunctionPairs;
189
190 // unique identifier per filter class
191 const char* _type;
192};
defines the cell type with cell-averaged quantities only (no linked cells).
Definition CouplingCell.h:29
const tarch::la::Vector< dim, double > & getMicroscopicMomentum() const
Definition CouplingCell.h:51
const tarch::la::Vector< dim, double > & getCurrentVelocity() const
Definition CouplingCell.h:92
const tarch::la::Vector< dim, double > & getMacroscopicMomentum() const
Definition CouplingCell.h:64
void setTemperature(const double &temperature)
Definition CouplingCell.h:96
void setMicroscopicMass(const double &mass)
Definition CouplingCell.h:42
const double & getPotentialEnergy() const
Definition CouplingCell.h:68
const double & getMicroscopicMass() const
Definition CouplingCell.h:45
void setMacroscopicMomentum(const tarch::la::Vector< dim, double > &momentum)
Definition CouplingCell.h:61
void setPotentialEnergy(const double &potentialEnergy)
Definition CouplingCell.h:71
const double & getTemperature() const
Definition CouplingCell.h:99
const double & getMacroscopicMass() const
Definition CouplingCell.h:58
void setMacroscopicMass(const double &mass)
Definition CouplingCell.h:55
void setCurrentVelocity(const tarch::la::Vector< dim, double > &velocity)
Definition CouplingCell.h:89
void setMicroscopicMomentum(const tarch::la::Vector< dim, double > &momentum)
Definition CouplingCell.h:48
Definition FilterInterface.h:32
std::vector< coupling::datastructures::CouplingCell< dim > * > _inputCells
Definition FilterInterface.h:181
Definition CellIndex.h:85
value_T get() const
Definition CellIndex.h:138
Definition Vector.h:24
Definition FilterPipeline.h:21
everything necessary for coupling operations, is defined in here
Definition AdditiveMomentumInsertion.h:15