OpenLB 1.7
Loading...
Searching...
No Matches
Public Member Functions | List of all members
olb::STLreader< T > Class Template Reference

#include <stlReader.h>

+ Inheritance diagram for olb::STLreader< T >:
+ Collaboration diagram for olb::STLreader< T >:

Public Member Functions

 STLreader (const std::string fName, T voxelSize, T stlSize=1, int method=2, bool verbose=false, T overlap=0., T max=0.)
 Constructs a new STLreader from a file.
 
 STLreader (const std::vector< std::vector< T > > meshPoints, T voxelSize, T stlSize=1, int method=2, bool verbose=false, T overlap=0., T max=0.)
 Constructs a new STLreader from a file.
 
 ~STLreader () override
 
bool operator() (bool output[], const T input[]) override
 Returns whether node is inside or not.
 
bool distance (T &distance, const Vector< T, 3 > &origin, const Vector< T, 3 > &direction, int iC=-1) override
 Computes distance to closest triangle intersection.
 
signedDistance (const Vector< T, 3 > &input) override
 Computes signed distance to closest triangle in direction of the surface normal.
 
Vector< T, 3 > surfaceNormal (const Vector< T, 3 > &pos, const T meshSize=0) override
 Finds and returns normal of the closest surface (triangle)
 
void print ()
 Prints console output.
 
void writeSTL (std::string stlName="")
 Writes STL mesh in Si units.
 
void writeOctree ()
 Writes Octree.
 
void setNormalsOutside ()
 Rearranges normals of triangles to point outside of geometry.
 
void setBoundaryInsideNodes ()
 Every octree leaf intersected by the STL will be part of the inside nodes.
 
Octree< T > * getTree () const
 Returns tree.
 
STLmesh< T > & getMesh ()
 Returns mesh.
 
- Public Member Functions inherited from olb::IndicatorF3D< T >
virtual Vector< T, 3 > & getMin ()
 
virtual Vector< T, 3 > & getMax ()
 
virtual bool distance (T &distance, const Vector< T, 3 > &origin, T precision, const Vector< T, 3 > &direction)
 
virtual bool distance (T &distance, const Vector< T, 3 > &origin, const Vector< T, 3 > &direction, T precision, T pitch)
 
virtual bool distance (T &distance, const Vector< T, 3 > &origin, const Vector< T, 3 > &direction, int iC=-1)
 
virtual bool distance (T &distance, const Vector< T, 3 > &origin)
 
virtual bool distance (T &distance, const T input[])
 
virtual bool normal (Vector< T, 3 > &normal, const Vector< T, 3 > &origin, const Vector< T, 3 > &direction, int iC=-1)
 returns true and the normal if there was one found for an given origin and direction
 
virtual bool rotOnAxis (Vector< T, 3 > &vec_rot, const Vector< T, 3 > &vec, const Vector< T, 3 > &axis, T &theta)
 Rotate vector around axis by angle theta.
 
virtual bool operator() (bool output[1], const T input[3])
 Returns true if input is inside the indicator.
 
virtual T signedDistance (const Vector< T, 3 > &input)
 Returns signed distance to the nearest point on the indicator surface.
 
virtual Vector< T, 3 > surfaceNormal (const Vector< T, 3 > &pos, const T meshSize)
 Return surface normal.
 
Vector< T, 3 > surfaceNormal (const Vector< T, 3 > &pos, const T meshSize, std::function< Vector< T, 3 >(const Vector< T, 3 > &)> transformPos)
 Return surface normal after possible translation and rotation.
 
bool isInsideBox (Vector< T, 3 > point)
 Returns true if point is inside a cube with corners _myMin and _myMax
 
virtual Vector< T, 3 > getSample (const std::function< T()> &randomness) const
 
- Public Member Functions inherited from olb::GenericF< T, S >
virtual ~GenericF ()=default
 
int getSourceDim () const
 read only access to member variable _m
 
int getTargetDim () const
 read only access to member variable _n
 
std::string & getName ()
 read and write access to name
 
std::string const & getName () const
 read only access to name
 
virtual bool operator() (T output[], const S input[])=0
 has to be implemented for 'every' derived class
 
bool operator() (T output[])
 wrapper that call the pure virtual operator() (T output[], const S input[]) from above
 
bool operator() (T output[], S input0)
 
bool operator() (T output[], S input0, S input1)
 
bool operator() (T output[], S input0, S input1, S input2)
 
bool operator() (T output[], S input0, S input1, S input2, S input3)
 

Additional Inherited Members

- Public Types inherited from olb::GenericF< T, S >
using targetType = T
 
using sourceType = S
 
- Public Attributes inherited from olb::GenericF< T, S >
std::shared_ptr< GenericF< T, S > > _ptrCalcC
 memory management, frees resouces (calcClass)
 
- Protected Member Functions inherited from olb::IndicatorF3D< T >
 IndicatorF3D ()
 
- Protected Member Functions inherited from olb::GenericF< T, S >
 GenericF (int targetDim, int sourceDim)
 
- Protected Attributes inherited from olb::IndicatorF3D< T >
Vector< T, 3 > _myMin
 
Vector< T, 3 > _myMax
 

Detailed Description

template<typename T>
class olb::STLreader< T >

Definition at line 214 of file stlReader.h.

Constructor & Destructor Documentation

◆ STLreader() [1/2]

template<typename T >
olb::STLreader< T >::STLreader ( const std::string fName,
T voxelSize,
T stlSize = 1,
int method = 2,
bool verbose = false,
T overlap = 0.,
T max = 0. )

Constructs a new STLreader from a file.

Parameters
fNameThe STL file name
voxelSizeVoxelsize in SI units
stlSizeConversion factor for STL (e.g. STL in mm stlSize=10^-3)
methodChoose indication method 0: fast, less stable 1: slow, more stable (for untight STLs)
verboseGet additional information.

Find center of tree and move by _voxelSize/4.

Create tree

Compute _myMin, _myMax such that they are the smallest (greatest) Voxel inside the STL.

Indicate nodes of the tree. (Inside/Outside)

Definition at line 713 of file stlReader.hh.

715 : _voxelSize(voxelSize),
716 _stlSize(stlSize),
717 _overlap(overlap),
718 _fName(fName),
719 _mesh(fName, stlSize),
720 _verbose(verbose),
721 clout(std::cout, "STLreader")
722{
723 this->getName() = "STLreader";
724
725 if (_verbose) {
726 clout << "Voxelizing ..." << std::endl;
727 }
728
729 Vector<T,3> extension = _mesh.getMax() - _mesh.getMin();
730 if ( util::nearZero(max) ) {
731 max = util::max(extension[0], util::max(extension[1], extension[2])) + _voxelSize;
732 }
733 int j = 0;
734 for (; _voxelSize * util::pow(2, j) < max; j++)
735 ;
736 Vector<T,3> center;
737 T radius = _voxelSize * util::pow(2, j - 1);
738
740 for (unsigned i = 0; i < 3; i++) {
741 center[i] = (_mesh.getMin()[i] + _mesh.getMax()[i]) / 2. - _voxelSize / 4.;
742 }
743
745 _tree = new Octree<T>(center, radius, &_mesh, j, _overlap);
746
748 for (int i = 0; i < 3; i++) {
749 this->_myMin[i] = center[i] + _voxelSize / 2.;
750 this->_myMax[i] = center[i] - _voxelSize / 2.;
751 }
752 for (int i = 0; i < 3; i++) {
753 while (this->_myMin[i] > _mesh.getMin()[i]) {
754 this->_myMin[i] -= _voxelSize;
755 }
756 while (this->_myMax[i] < _mesh.getMax()[i]) {
757 this->_myMax[i] += _voxelSize;
758 }
759 this->_myMax[i] -= _voxelSize;
760 this->_myMin[i] += _voxelSize;
761 }
762
764 switch (method) {
765 case 1:
766 indicate1();
767 break;
768 case 3:
769 indicate3();
770 break;
771 case 4:
772 indicate2_Xray();
773 break;
774 case 5:
775 indicate2_Yray();
776 break;
777 default:
778 indicate2();
779 break;
780 }
781
782 if (_verbose) {
783 print();
784 }
785 if (_verbose) {
786 clout << "Voxelizing ... OK" << std::endl;
787 }
788}
std::string & getName()
read and write access to name
Definition genericF.hh:51
void print()
Prints console output.
Pack< T > max(Pack< T > rhs, Pack< T > lhs)
Definition 256.h:416
cpu::simd::Pack< T > max(cpu::simd::Pack< T > rhs, cpu::simd::Pack< T > lhs)
Definition pack.h:130
cpu::simd::Pack< T > pow(cpu::simd::Pack< T > base, cpu::simd::Pack< T > exp)
Definition pack.h:112
bool nearZero(const ADf< T, DIM > &a)
Definition aDiff.h:1087

References olb::IndicatorF3D< T >::_myMax, olb::IndicatorF3D< T >::_myMin, olb::GenericF< T, S >::getName(), olb::util::max(), olb::util::nearZero(), olb::util::pow(), and olb::STLreader< T >::print().

+ Here is the call graph for this function:

◆ STLreader() [2/2]

template<typename T >
olb::STLreader< T >::STLreader ( const std::vector< std::vector< T > > meshPoints,
T voxelSize,
T stlSize = 1,
int method = 2,
bool verbose = false,
T overlap = 0.,
T max = 0. )

Constructs a new STLreader from a file.

Parameters
fNameThe STL file name
voxelSizeVoxelsize in SI units
stlSizeConversion factor for STL (e.g. STL in mm stlSize=10^-3)
methodChoose indication method 0: fast, less stable 1: slow, more stable (for untight STLs)
verboseGet additional information.

Find center of tree and move by _voxelSize/4.

Create tree

Compute _myMin, _myMax such that they are the smallest (greatest) Voxel inside the STL.

Definition at line 794 of file stlReader.hh.

796 : _voxelSize(voxelSize),
797 _stlSize(stlSize),
798 _overlap(overlap),
799 _fName("meshPoints.stl"),
800 _mesh(meshPoints, stlSize),
801 _verbose(verbose),
802 clout(std::cout, "STLreader")
803{
804 this->getName() = "STLreader";
805
806 if (_verbose) {
807 clout << "Voxelizing ..." << std::endl;
808 }
809
810 Vector<T,3> extension = _mesh.getMax() - _mesh.getMin();
811 if ( util::nearZero(max) ) {
812 max = util::max(extension[0], util::max(extension[1], extension[2])) + _voxelSize;
813 }
814 int j = 0;
815 for (; _voxelSize * util::pow(2, j) < max; j++)
816 ;
817 Vector<T,3> center;
818 T radius = _voxelSize * util::pow(2, j - 1);
819
821 for (unsigned i = 0; i < 3; i++) {
822 center[i] = (_mesh.getMin()[i] + _mesh.getMax()[i]) / 2. - _voxelSize / 4.;
823 }
824
826
827 _tree = new Octree<T>(center, radius, &_mesh, j, _overlap);
828
830 for (int i = 0; i < 3; i++) {
831 this->_myMin[i] = center[i] + _voxelSize / 2.;
832 this->_myMax[i] = center[i] - _voxelSize / 2.;
833 }
834 for (int i = 0; i < 3; i++) {
835 while (this->_myMin[i] > _mesh.getMin()[i]) {
836 this->_myMin[i] -= _voxelSize;
837 }
838 while (this->_myMax[i] < _mesh.getMax()[i]) {
839 this->_myMax[i] += _voxelSize;
840 }
841 this->_myMax[i] -= _voxelSize;
842 this->_myMin[i] += _voxelSize;
843 }
844 //automaticly choose the method with minimum extension in its direction
845
846 /*if(extension[0] == std::min_element(extension.begin(), extension.end())){
847 method = 4;
848 }
849 else if(extension[1] == std::min_element(extension.begin(), extension.end())){
850 method = 5;
851 }
852 else if(extension[2] == std::min_element(extension.begin(), extension.end())){
853 method = 0;
854 }
855 */
856
857
858 // Indicate nodes of the tree. (Inside/Outside)
859 switch (method) {
860 case 1:
861 indicate1();
862 break;
863 case 3:
864 indicate3();
865 break;
866 case 4:
867 indicate2_Xray();
868 break;
869 case 5:
870 indicate2_Yray();
871 break;
872 default:
873 indicate2();
874 break;
875 }
876
878
879 if (_verbose) {
880 print();
881 }
882 if (_verbose) {
883 clout << "Voxelizing ... OK" << std::endl;
884 }
885}
void setNormalsOutside()
Rearranges normals of triangles to point outside of geometry.

References olb::IndicatorF3D< T >::_myMax, olb::IndicatorF3D< T >::_myMin, olb::GenericF< T, S >::getName(), olb::util::max(), olb::util::nearZero(), olb::util::pow(), olb::STLreader< T >::print(), and olb::STLreader< T >::setNormalsOutside().

+ Here is the call graph for this function:

◆ ~STLreader()

template<typename T >
olb::STLreader< T >::~STLreader ( )
override

Definition at line 888 of file stlReader.hh.

889{
890 delete _tree;
891}

Member Function Documentation

◆ distance()

template<typename T >
bool olb::STLreader< T >::distance ( T & distance,
const Vector< T, 3 > & origin,
const Vector< T, 3 > & direction,
int iC = -1 )
override

Computes distance to closest triangle intersection.

Definition at line 1241 of file stlReader.hh.

1243{
1244 Octree<T>* node = nullptr;
1245 Vector<T,3> dir(direction);
1246 dir = normalize(dir);
1247 Vector<T,3> extends = _mesh.getMax() - _mesh.getMin();
1248 Vector<T,3> pt(origin);
1249 Vector<T,3> q;
1250 Vector<T,3> s;
1251 Vector<T,3> center = _mesh.getMin() + 1 / 2. * extends;
1252 T step = _voxelSize / 1000., a = 0;
1253
1254 for (int i = 0; i < 3; i++) {
1255 extends[i] /= 2.;
1256 }
1257
1258 if (!(_mesh.getMin()[0] < origin[0] && origin[0] < _mesh.getMax()[0]
1259 && _mesh.getMin()[1] < origin[1] && origin[1] < _mesh.getMax()[1]
1260 && _mesh.getMin()[2] < origin[2] && origin[2] < _mesh.getMax()[2])) {
1261 T t = T(), d = T();
1262 bool foundQ = false;
1263
1264 if (dir[0] > 0) {
1265 d = _mesh.getMin()[0];
1266 t = (d - origin[0]) / dir[0];
1267 pt[0] = origin[0] + (t + step) * dir[0];
1268 pt[1] = origin[1] + (t + step) * dir[1];
1269 pt[2] = origin[2] + (t + step) * dir[2];
1270
1271 if (_mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]
1272 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1273 foundQ = true;
1274 }
1275 }
1276 else if (dir[0] < 0) {
1277 d = _mesh.getMax()[0];
1278 t = (d - origin[0]) / dir[0];
1279 pt[0] = origin[0] + (t + step) * dir[0];
1280 pt[1] = origin[1] + (t + step) * dir[1];
1281 pt[2] = origin[2] + (t + step) * dir[2];
1282 if (_mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]
1283 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1284 foundQ = true;
1285 }
1286 }
1287
1288 if (dir[1] > 0 && !foundQ) {
1289 d = _mesh.getMin()[1];
1290 t = (d - origin[1]) / dir[1];
1291 pt[0] = origin[0] + (t + step) * dir[0];
1292 pt[1] = origin[1] + (t + step) * dir[1];
1293 pt[2] = origin[2] + (t + step) * dir[2];
1294 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1295 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1296 foundQ = true;
1297 }
1298 }
1299 else if (dir[1] < 0 && !foundQ) {
1300 d = _mesh.getMax()[1];
1301 t = (d - origin[1]) / dir[1];
1302 pt[0] = origin[0] + (t + step) * dir[0];
1303 pt[1] = origin[1] + (t + step) * dir[1];
1304 pt[2] = origin[2] + (t + step) * dir[2];
1305 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1306 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1307 foundQ = true;
1308 }
1309 }
1310
1311 if (dir[2] > 0 && !foundQ) {
1312 d = _mesh.getMin()[2];
1313 t = (d - origin[2]) / dir[2];
1314 pt[0] = origin[0] + (t + step) * dir[0];
1315 pt[1] = origin[1] + (t + step) * dir[1];
1316 pt[2] = origin[2] + (t + step) * dir[2];
1317 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1318 && _mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]) {
1319 foundQ = true;
1320 }
1321 }
1322 else if (dir[2] < 0 && !foundQ) {
1323 d = _mesh.getMax()[2];
1324 t = (d - origin[2]) / dir[2];
1325 pt[0] = origin[0] + (t + step) * dir[0];
1326 pt[1] = origin[1] + (t + step) * dir[1];
1327 pt[2] = origin[2] + (t + step) * dir[2];
1328 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1329 && _mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]) {
1330 foundQ = true;
1331 }
1332 }
1333
1334 if (!foundQ) {
1335 return false;
1336 }
1337 }
1338
1339 while ((util::fabs(pt[0] - center[0]) < extends[0])
1340 && (util::fabs(pt[1] - center[1]) < extends[1])
1341 && (util::fabs(pt[2] - center[2]) < extends[2])) {
1342 node = _tree->find(pt);
1343 if (node->closestIntersection(Vector<T,3>(origin), dir, q, a)) {
1344 Vector<T,3> vek(q - Vector<T,3>(origin));
1345 distance = norm(vek);
1346 return true;
1347 }
1348 else {
1349 Octree<T>* tmpNode = _tree->find(pt);
1350 tmpNode->intersectRayNode(pt, dir, s);
1351 for (int i = 0; i < 3; i++) {
1352 pt[i] = s[i] + step * dir[i];
1353 }
1354 }
1355 }
1356
1357 return false;
1358}
bool distance(T &distance, const Vector< T, 3 > &origin, const Vector< T, 3 > &direction, int iC=-1) override
Computes distance to closest triangle intersection.
platform_constant Fraction t[Q]
platform_constant Fraction s[Q]
constexpr int q() any_platform
constexpr int d() any_platform
cpu::simd::Pack< T > fabs(cpu::simd::Pack< T > value)
Definition pack.h:106
constexpr T norm(const ScalarVector< T, D, IMPL > &a)
Euclidean vector norm.
constexpr Vector< T, D > normalize(const ScalarVector< T, D, IMPL > &a, T scale=T{1})
Definition vector.h:245

References olb::Octree< T >::closestIntersection(), olb::util::fabs(), olb::Octree< T >::find(), olb::Octree< T >::intersectRayNode(), olb::norm(), and olb::normalize().

+ Here is the call graph for this function:

◆ getMesh()

template<typename T >
STLmesh< T > & olb::STLreader< T >::getMesh ( )
inline

Returns mesh.

Definition at line 344 of file stlReader.h.

345 {
346 return _mesh;
347 };

◆ getTree()

template<typename T >
Octree< T > * olb::STLreader< T >::getTree ( ) const
inline

Returns tree.

Definition at line 337 of file stlReader.h.

338 {
339 return _tree;
340 };

◆ operator()()

template<typename T >
bool olb::STLreader< T >::operator() ( bool output[],
const T input[] )
override

Returns whether node is inside or not.

Definition at line 1227 of file stlReader.hh.

1228{
1229 output[0] = false;
1230 T coords = _tree->getRadius();
1231 Vector<T,3> c(_tree->getCenter());
1232 if (c[0] - coords < input[0] && input[0] < c[0] + coords && c[1] - coords < input[1]
1233 && input[1] < c[1] + coords && c[2] - coords < input[2] && input[2] < c[2] + coords) {
1234 std::vector<T> tmp(input, input + 3);
1235 output[0] = _tree->find(tmp)->getInside();
1236 }
1237 return true;
1238}
platform_constant int c[Q][D]

◆ print()

template<typename T >
void olb::STLreader< T >::print ( )

Prints console output.

Definition at line 1494 of file stlReader.hh.

1495{
1496 _mesh.print();
1497 _tree->print();
1498 clout << "voxelSize=" << _voxelSize << "; stlSize=" << _stlSize << std::endl;
1499 clout << "minPhysR(VoxelMesh)=(" << this->_myMin[0] << "," << this->_myMin[1]
1500 << "," << this->_myMin[2] << ")";
1501 clout << "; maxPhysR(VoxelMesh)=(" << this->_myMax[0] << ","
1502 << this->_myMax[1] << "," << this->_myMax[2] << ")" << std::endl;
1503}
+ Here is the caller graph for this function:

◆ setBoundaryInsideNodes()

template<typename T >
void olb::STLreader< T >::setBoundaryInsideNodes ( )

Every octree leaf intersected by the STL will be part of the inside nodes.

Artificially enlarges all details that would otherwise be cut off by the voxelSize.

Definition at line 1545 of file stlReader.hh.

1546{
1547 std::vector<Octree<T>*> leafs;
1548 _tree->getLeafs(leafs);
1549 for (auto it = leafs.begin(); it != leafs.end(); ++it) {
1550 if ((*it)->getBoundaryNode()) {
1551 (*it)->setInside(true);
1552 }
1553 }
1554}

◆ setNormalsOutside()

template<typename T >
void olb::STLreader< T >::setNormalsOutside ( )

Rearranges normals of triangles to point outside of geometry.

Definition at line 1523 of file stlReader.hh.

1524{
1525 unsigned int noTris = _mesh.triangleSize();
1526 Vector<T,3> center;
1527 //Octree<T>* node = nullptr;
1528 for (unsigned int i = 0; i < noTris; i++) {
1529 center = _mesh.getTri(i).getCenter();
1530 if (_tree->find(
1531 center + _mesh.getTri(i).normal * util::sqrt(3.) * _voxelSize)->getInside()) {
1532 // cout << "Wrong direction" << std::endl;
1533 Vector<T,3> pt(_mesh.getTri(i).point[0].coords);
1534 _mesh.getTri(i).point[0].coords = _mesh.getTri(i).point[2].coords;
1535 _mesh.getTri(i).point[2].coords = pt;
1536 _mesh.getTri(i).init();
1537 // _mesh.getTri(i).getNormal()[0] *= -1.;
1538 // _mesh.getTri(i).getNormal()[1] *= -1.;
1539 // _mesh.getTri(i).getNormal()[2] *= -1.;
1540 }
1541 }
1542}
cpu::simd::Pack< T > sqrt(cpu::simd::Pack< T > value)
Definition pack.h:100

References olb::util::sqrt().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ signedDistance()

template<typename T >
T olb::STLreader< T >::signedDistance ( const Vector< T, 3 > & input)
override

Computes signed distance to closest triangle in direction of the surface normal.

Definition at line 1467 of file stlReader.hh.

1468{
1469 Vector<T,3> distance;
1470 Vector<T,3> closestPointOnSurface;
1471 T distanceNorm = std::numeric_limits<T>::max();
1472
1473 for (const STLtriangle<T>& triangle : _mesh.getTriangles()) {
1474 const PhysR<T,3> currPointOnTriangle = triangle.closestPtPointTriangle(input);
1475 const Vector<T,3> currDistance = input - currPointOnTriangle;
1476 T currDistanceNorm = norm_squared(currDistance);
1477 if (distanceNorm > currDistanceNorm) {
1478 distanceNorm = currDistanceNorm;
1479 distance = currDistance;
1480 closestPointOnSurface = currPointOnTriangle;
1481 }
1482 }
1483
1484 return util::sqrt(distanceNorm) * evalSignForSignedDistance(input);
1485}
T triangle(Vector< T, 2 > p, Vector< T, 2 > a, Vector< T, 2 > b, Vector< T, 2 > c) any_platform
Exact signed distance to the surface of two-dimensional triangle.
Definition sdf.h:143
constexpr T norm_squared(const ScalarVector< T, D, IMPL > &a)
Squared euclidean vector norm.

References olb::norm_squared(), and olb::util::sqrt().

+ Here is the call graph for this function:

◆ surfaceNormal()

template<typename T >
Vector< T, 3 > olb::STLreader< T >::surfaceNormal ( const Vector< T, 3 > & pos,
const T meshSize = 0 )
override

Finds and returns normal of the closest surface (triangle)

Definition at line 1488 of file stlReader.hh.

1489{
1490 return evalSurfaceNormal(pos);
1491}

◆ writeOctree()

template<typename T >
void olb::STLreader< T >::writeOctree ( )

Writes Octree.

Definition at line 1506 of file stlReader.hh.

1507{
1508 _tree->write(_fName);
1509}

◆ writeSTL()

template<typename T >
void olb::STLreader< T >::writeSTL ( std::string stlName = "")

Writes STL mesh in Si units.

Definition at line 1512 of file stlReader.hh.

1513{
1514 if (stlName == "") {
1515 _mesh.write(_fName);
1516 }
1517 else {
1518 _mesh.write(stlName);
1519 }
1520}

The documentation for this class was generated from the following files: