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
7namespace viennals {
8
9using namespace viennacore;
10
14template <class T> class Extrude {
15 SmartPointer<Domain<T, 2>> inputLevelSet = nullptr;
16 SmartPointer<Domain<T, 3>> outputLevelSet = nullptr;
17 Vec2D<T> extent = {0., 0.};
18 int extrudeDim = 0;
19 std::array<BoundaryConditionEnum<3>, 3> boundaryConds;
20
21public:
23 Extrude(SmartPointer<Domain<T, 2>> passedInputLS,
24 SmartPointer<Domain<T, 3>> passedOutputLS, Vec2D<T> passedExtent,
25 const int passedExtrudeDim,
26 std::array<BoundaryConditionEnum<3>, 3> passedBoundaryConds)
27 : inputLevelSet(passedInputLS), outputLevelSet(passedOutputLS),
28 extent(passedExtent), extrudeDim(passedExtrudeDim),
29 boundaryConds(passedBoundaryConds) {}
30
31 void setInputLevelSet(SmartPointer<Domain<T, 2>> passedInputLS) {
32 inputLevelSet = passedInputLS;
33 }
34
35 // The 3D output LS will be overwritten by the extruded LS
36 void setOutputLevelSet(SmartPointer<Domain<T, 3>> &passedOutputLS) {
37 outputLevelSet = passedOutputLS;
38 }
39
40 // Set the min and max extent in the extruded dimension
41 void setExtent(Vec2D<T> passedExtent) { extent = passedExtent; }
42
43 // Set which index of the added dimension (x: 0, y: 1, z: 2)
44 void setExtrudeDimension(const int passedExtrudeDim) {
45 extrudeDim = passedExtrudeDim;
46 }
47
49 std::array<BoundaryConditionEnum<3>, 3> passedBoundaryConds) {
50 boundaryConds = passedBoundaryConds;
51 }
52
53 void setBoundaryConditions(BoundaryConditionEnum<3> passedBoundaryConds[3]) {
54 for (int i = 0; i < 3; i++)
55 boundaryConds[i] = passedBoundaryConds[i];
56 }
57
58 void apply() {
59 if (inputLevelSet == nullptr) {
60 Logger::getInstance()
61 .addWarning("No input Level Set supplied to Extrude! Not converting.")
62 .print();
63 }
64 if (outputLevelSet == nullptr) {
65 Logger::getInstance()
66 .addWarning(
67 "No output Level Set supplied to Extrude! Not converting.")
68 .print();
69 return;
70 }
71
72 // x and y of the input LS get transformed to these indices
73 const auto extrudeDims = getExtrudeDims();
74
75 // create new domain based on 2D extent
76 {
77 const T gridDelta = inputLevelSet->getGrid().getGridDelta();
78 auto minBounds = inputLevelSet->getGrid().getMinBounds();
79 auto maxBounds = inputLevelSet->getGrid().getMaxBounds();
80
81 double domainBounds[2 * 3];
82 domainBounds[2 * extrudeDim] = extent[0];
83 domainBounds[2 * extrudeDim + 1] = extent[1];
84 domainBounds[2 * extrudeDims[0]] = gridDelta * minBounds[0];
85 domainBounds[2 * extrudeDims[0] + 1] = gridDelta * maxBounds[0];
86 domainBounds[2 * extrudeDims[1]] = gridDelta * minBounds[1];
87 domainBounds[2 * extrudeDims[1] + 1] = gridDelta * maxBounds[1];
88
89 auto tmpLevelSet = SmartPointer<Domain<T, 3>>::New(
90 domainBounds, boundaryConds.data(), gridDelta);
91 outputLevelSet->deepCopy(tmpLevelSet);
92 }
93
94 auto surface = SmartPointer<Mesh<T>>::New();
95 ToSurfaceMesh<T, 2>(inputLevelSet, surface).apply();
96
97 auto &lines = surface->template getElements<2>();
98 auto &nodes = surface->getNodes();
99 const unsigned numNodes = nodes.size();
100
101 // add new nodes shifted by the extent
102 for (unsigned i = 0; i < numNodes; i++) {
103 nodes[i][extrudeDims[1]] = nodes[i][1];
104 nodes[i][extrudeDims[0]] = nodes[i][0];
105 nodes[i][extrudeDim] = extent[1];
106
107 nodes.push_back(nodes[i]);
108
109 nodes[i][extrudeDim] = extent[0];
110 }
111
112 // add triangles in places of lines
113 for (unsigned i = 0; i < lines.size(); i++) {
114 std::array<unsigned, 3> triangle = {lines[i][1], lines[i][0],
115 lines[i][0] + numNodes};
116 if (extrudeDim == 1) {
117 std::swap(triangle[0], triangle[2]);
118 }
119 surface->insertNextTriangle(triangle);
120 triangle[0] = lines[i][0] + numNodes;
121 triangle[1] = lines[i][1] + numNodes;
122 triangle[2] = lines[i][1];
123 if (extrudeDim == 1) {
124 std::swap(triangle[0], triangle[2]);
125 }
126 surface->insertNextTriangle(triangle);
127 }
128 surface->template getElements<2>().clear(); // remove lines
129
130 FromSurfaceMesh<T, 3>(outputLevelSet, surface).apply();
131 }
132
133private:
134 inline std::array<unsigned, 2> getExtrudeDims() const {
135 assert(extrudeDim < 3);
136 if (extrudeDim == 0) {
137 return {1, 2};
138 } else if (extrudeDim == 1) {
139 return {0, 2};
140 } else {
141 return {0, 1};
142 }
143 }
144};
145
146} // namespace viennals
Class containing all information about the level set, including the dimensions of the domain,...
Definition lsDomain.hpp:28
Extrudes a 2D Level Set into a 3D domain. The axis in which should be extruded can be set and boundar...
Definition lsExtrude.hpp:14
void setExtent(Vec2D< T > passedExtent)
Definition lsExtrude.hpp:41
Extrude()
Definition lsExtrude.hpp:22
Extrude(SmartPointer< Domain< T, 2 > > passedInputLS, SmartPointer< Domain< T, 3 > > passedOutputLS, Vec2D< T > passedExtent, const int passedExtrudeDim, std::array< BoundaryConditionEnum< 3 >, 3 > passedBoundaryConds)
Definition lsExtrude.hpp:23
void setBoundaryConditions(BoundaryConditionEnum< 3 > passedBoundaryConds[3])
Definition lsExtrude.hpp:53
void apply()
Definition lsExtrude.hpp:58
void setExtrudeDimension(const int passedExtrudeDim)
Definition lsExtrude.hpp:44
void setOutputLevelSet(SmartPointer< Domain< T, 3 > > &passedOutputLS)
Definition lsExtrude.hpp:36
void setInputLevelSet(SmartPointer< Domain< T, 2 > > passedInputLS)
Definition lsExtrude.hpp:31
void setBoundaryConditions(std::array< BoundaryConditionEnum< 3 >, 3 > passedBoundaryConds)
Definition lsExtrude.hpp:48
Construct a level set from an explicit mesh.
Definition lsFromSurfaceMesh.hpp:15
void apply()
Definition lsFromSurfaceMesh.hpp:250
Extract an explicit Mesh<> instance from an lsDomain. The interface is then described by explciit sur...
Definition lsToSurfaceMesh.hpp:20
void apply()
Definition lsToSurfaceMesh.hpp:44
Definition lsAdvect.hpp:46
typename hrleGrid< D >::boundaryType BoundaryConditionEnum
Definition lsDomain.hpp:24
double T
Definition pyWrap.cpp:63