ViennaLS
Loading...
Searching...
No Matches
lsExtrude.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <lsDomain.hpp>
5#include <lsToSurfaceMesh.hpp>
6
7#include <hrleSparseIterator.hpp>
8
9namespace viennals {
10
11using namespace viennacore;
12
16template <class T> class Extrude {
17 using hrleIndexType = viennahrle::IndexType;
18 SmartPointer<Domain<T, 2>> inputLevelSet = nullptr;
19 SmartPointer<Domain<T, 3>> outputLevelSet = nullptr;
20 Vec2D<T> extent = {0., 0.};
21 int extrudeDim = 0;
22 std::array<BoundaryConditionEnum, 3> boundaryConds = {};
23
24public:
25 Extrude() = default;
26
27 Extrude(SmartPointer<Domain<T, 2>> passedInputLS,
28 SmartPointer<Domain<T, 3>> passedOutputLS, Vec2D<T> passedExtent,
29 const int passedExtrudeDim = 2,
30 BoundaryConditionEnum passedBoundaryCond =
31 BoundaryConditionEnum::INFINITE_BOUNDARY)
32 : inputLevelSet(passedInputLS), outputLevelSet(passedOutputLS),
33 extent(passedExtent), extrudeDim(passedExtrudeDim),
34 boundaryConds{{passedInputLS->getGrid().getBoundaryConditions(0),
35 passedInputLS->getGrid().getBoundaryConditions(1),
36 passedBoundaryCond}} {}
37
38 Extrude(SmartPointer<Domain<T, 2>> passedInputLS,
39 SmartPointer<Domain<T, 3>> passedOutputLS, Vec2D<T> passedExtent,
40 const int passedExtrudeDim,
41 std::array<BoundaryConditionEnum, 3> passedBoundaryConds)
42 : inputLevelSet(passedInputLS), outputLevelSet(passedOutputLS),
43 extent(passedExtent), extrudeDim(passedExtrudeDim),
44 boundaryConds(passedBoundaryConds) {}
45
46 void setInputLevelSet(SmartPointer<Domain<T, 2>> passedInputLS) {
47 inputLevelSet = passedInputLS;
48 }
49
50 // The 3D output LS will be overwritten by the extruded LS
51 void setOutputLevelSet(SmartPointer<Domain<T, 3>> &passedOutputLS) {
52 outputLevelSet = passedOutputLS;
53 }
54
55 // Set the min and max extent in the extruded dimension
56 void setExtent(Vec2D<T> passedExtent) { extent = passedExtent; }
57
58 // Set which index of the added dimension (x: 0, y: 1, z: 2)
59 void setExtrudeDimension(const int passedExtrudeDim) {
60 extrudeDim = passedExtrudeDim;
61 }
62
64 std::array<BoundaryConditionEnum, 3> passedBoundaryConds) {
65 boundaryConds = passedBoundaryConds;
66 }
67
68 void setBoundaryConditions(BoundaryConditionEnum passedBoundaryConds[3]) {
69 for (int i = 0; i < 3; i++)
70 boundaryConds[i] = passedBoundaryConds[i];
71 }
72
73 void apply() {
74 if (inputLevelSet == nullptr) {
75 Logger::getInstance()
76 .addError("No input Level Set supplied to Extrude.")
77 .print();
78 }
79 if (outputLevelSet == nullptr) {
80 Logger::getInstance()
81 .addError("No output Level Set supplied to Extrude.")
82 .print();
83 return;
84 }
85
86 // x and y of the input LS get transformed to these indices
87 const auto extrudeDims = getExtrudeDims();
88
89 std::vector<std::pair<viennahrle::Index<3>, T>> points3D;
90
91 auto &domain2D = inputLevelSet->getDomain();
92 auto &grid2D = inputLevelSet->getGrid();
93 const T gridDelta = grid2D.getGridDelta();
94 auto minBounds = grid2D.getMinBounds();
95 auto maxBounds = grid2D.getMaxBounds();
96
97 double domainBounds[2 * 3];
98 domainBounds[2 * extrudeDim] = extent[0];
99 domainBounds[2 * extrudeDim + 1] = extent[1];
100 domainBounds[2 * extrudeDims[0]] = gridDelta * minBounds[0];
101 domainBounds[2 * extrudeDims[0] + 1] = gridDelta * maxBounds[0];
102 domainBounds[2 * extrudeDims[1]] = gridDelta * minBounds[1];
103 domainBounds[2 * extrudeDims[1] + 1] = gridDelta * maxBounds[1];
104
105 auto tmpLevelSet = SmartPointer<Domain<T, 3>>::New(
106 domainBounds, boundaryConds.data(), gridDelta);
107 outputLevelSet->deepCopy(tmpLevelSet);
108
109 for (viennahrle::SparseIterator<typename Domain<T, 2>::DomainType> it(
110 domain2D);
111 !it.isFinished(); ++it) {
112 if (!it.isDefined())
113 continue;
114
115 const auto index2D = it.getStartIndices();
116 const T value = it.getValue();
117
118 const hrleIndexType zStart = std::floor(extent[0] / gridDelta) - 1;
119 const hrleIndexType zEnd = std::ceil(extent[1] / gridDelta) + 1;
120 for (hrleIndexType z = zStart; z <= zEnd; ++z) {
121 viennahrle::Index<3> index3D;
122 index3D[0] = index2D[0];
123 index3D[1] = index2D[1];
124 index3D[2] = z;
125
126 points3D.emplace_back(index3D, value);
127 }
128 }
129
130 outputLevelSet->insertPoints(points3D);
131 outputLevelSet->finalize(2);
132 }
133
134private:
135 inline std::array<unsigned, 2> getExtrudeDims() const {
136 assert(extrudeDim < 3);
137 if (extrudeDim == 0) {
138 return {1, 2};
139 } else if (extrudeDim == 1) {
140 return {0, 2};
141 } else {
142 return {0, 1};
143 }
144 }
145};
146
147} // namespace viennals
double T
Definition Epitaxy.cpp:10
Class containing all information about the level set, including the dimensions of the domain,...
Definition lsDomain.hpp:27
viennahrle::Domain< T, D > DomainType
Definition lsDomain.hpp:32
void setBoundaryConditions(BoundaryConditionEnum passedBoundaryConds[3])
Definition lsExtrude.hpp:68
void setExtent(Vec2D< T > passedExtent)
Definition lsExtrude.hpp:56
void setBoundaryConditions(std::array< BoundaryConditionEnum, 3 > passedBoundaryConds)
Definition lsExtrude.hpp:63
Extrude(SmartPointer< Domain< T, 2 > > passedInputLS, SmartPointer< Domain< T, 3 > > passedOutputLS, Vec2D< T > passedExtent, const int passedExtrudeDim, std::array< BoundaryConditionEnum, 3 > passedBoundaryConds)
Definition lsExtrude.hpp:38
void apply()
Definition lsExtrude.hpp:73
Extrude(SmartPointer< Domain< T, 2 > > passedInputLS, SmartPointer< Domain< T, 3 > > passedOutputLS, Vec2D< T > passedExtent, const int passedExtrudeDim=2, BoundaryConditionEnum passedBoundaryCond=BoundaryConditionEnum::INFINITE_BOUNDARY)
Definition lsExtrude.hpp:27
void setExtrudeDimension(const int passedExtrudeDim)
Definition lsExtrude.hpp:59
void setOutputLevelSet(SmartPointer< Domain< T, 3 > > &passedOutputLS)
Definition lsExtrude.hpp:51
void setInputLevelSet(SmartPointer< Domain< T, 2 > > passedInputLS)
Definition lsExtrude.hpp:46
Definition lsAdvect.hpp:36
viennahrle::BoundaryType BoundaryConditionEnum
Definition lsDomain.hpp:23