preCICE v3.1.2
Loading...
Searching...
No Matches
RadialGeoMultiscaleMapping.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <vector>
4#include "logging/Logger.hpp"
5#include "mapping/Mapping.hpp"
6
7namespace precice {
8namespace mapping {
9
12public:
19 enum struct MultiscaleType {
20 SPREAD,
22 };
23 enum struct MultiscaleAxis {
24 X = 0,
25 Y = 1,
26 Z = 2
27 };
28
37 RadialGeoMultiscaleMapping(Constraint constraint, int dimensions, MultiscaleType type, MultiscaleAxis axis);
38
40 void computeMapping() override;
41
43 void clear() override;
44
45 void tagMeshFirstRound() override;
46 void tagMeshSecondRound() override;
47
49 std::string getName() const final override;
50
51protected:
53 void mapConservative(const time::Sample &inData, Eigen::VectorXd &outData) override;
54
56 void mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData) override;
57
58private:
59 mutable logging::Logger _log{"mapping::RadialGeoMultiscaleMapping"};
60
63
66
70
73};
74
75} // namespace mapping
76} // namespace precice
Abstract base class for mapping of data from one mesh to another.
Definition Mapping.hpp:15
Constraint
Specifies additional constraints for a mapping.
Definition Mapping.hpp:29
Geometric multiscale mapping in radial direction.
MultiscaleType
Geometric multiscale type of the mapping (spread or collect).
RadialGeoMultiscaleMapping(Constraint constraint, int dimensions, MultiscaleType type, MultiscaleAxis axis)
Constructor.
std::vector< size_t > _vertexCounter
counts number of vertices between midpoints for averaging
MultiscaleAxis _axis
main axis along which radial geometric multiscale coupling happens
void tagMeshFirstRound() override
Method used by partition. Tags vertices that could be owned by this rank.
std::vector< size_t > _vertexIndicesSpread
computed vertex indices to map data from input vertices to output vertices and vice versa
std::string getName() const final override
Returns name of the mapping.
MultiscaleType _type
type of mapping, namely spread or collect
void clear() override
Removes a computed mapping.
void tagMeshSecondRound() override
Method used by partition. Tags vertices that can be filtered out.
void mapConservative(const time::Sample &inData, Eigen::VectorXd &outData) override
Maps data using a conservative constraint.
void computeMapping() override
Takes care of compute-heavy operations needed only once to set up the mapping.
void mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData) override
Maps data using a consistent constraint.
Main namespace of the precice library.