Skip to content
Snippets Groups Projects
Commit eae8a654 authored by Christian Godenschwager's avatar Christian Godenschwager
Browse files

Add simple StructuredBlockForest initialization functions

parent 6985b1e0
No related merge requests found
...@@ -21,6 +21,10 @@ ...@@ -21,6 +21,10 @@
#pragma once #pragma once
#include "BlockExclusion.h"
#include "mesh/MeshOperations.h"
#include "mesh/distance_octree/DistanceOctree.h"
#include "blockforest/StructuredBlockForest.h" #include "blockforest/StructuredBlockForest.h"
#include "blockforest/SetupBlockForest.h" #include "blockforest/SetupBlockForest.h"
#include "blockforest/loadbalancing/StaticCurve.h" #include "blockforest/loadbalancing/StaticCurve.h"
...@@ -110,5 +114,37 @@ private: ...@@ -110,5 +114,37 @@ private:
}; };
template< typename MeshType >
shared_ptr<StructuredBlockForest> createStructuredBlockStorageInsideMesh( const shared_ptr< mesh::DistanceOctree< MeshType > > & distanceOctree, const real_t dx, const Vector3<uint_t> & blockSize )
{
ComplexGeometryStructuredBlockforestCreator creator( distanceOctree->getAABB(), Vector3<real_t>(dx), makeExcludeMeshExterior( distanceOctree, dx ) );
return creator.createStructuredBlockForest( blockSize );
}
template< typename MeshType >
shared_ptr<StructuredBlockForest> createStructuredBlockStorageOutsideMesh( const AABB & aabb, const shared_ptr< mesh::DistanceOctree< MeshType > > & distanceOctree, const real_t dx, const Vector3<uint_t> & blockSize )
{
ComplexGeometryStructuredBlockforestCreator creator( aabb, Vector3<real_t>(dx), makeExcludeMeshInterior( distanceOctree, dx ) );
return creator.createStructuredBlockForest( blockSize );
}
template< typename MeshType >
shared_ptr<StructuredBlockForest> createStructuredBlockStorageInsideMesh( const shared_ptr< mesh::DistanceOctree< MeshType > > & distanceOctree, const real_t dx, const uint_t targetNumRootBlocks )
{
ComplexGeometryStructuredBlockforestCreator creator( distanceOctree->getAABB(), Vector3<real_t>(dx), makeExcludeMeshExterior( distanceOctree, dx ) );
return creator.createStructuredBlockForest( targetNumRootBlocks );
}
template< typename MeshType >
shared_ptr<StructuredBlockForest> createStructuredBlockStorageOutsideMesh( const AABB & aabb, const shared_ptr< mesh::DistanceOctree< MeshType > > & distanceOctree, const real_t dx, const uint_t targetNumRootBlocks )
{
ComplexGeometryStructuredBlockforestCreator creator( aabb, Vector3<real_t>(dx), makeExcludeMeshInterior( distanceOctree, dx ) );
return creator.createStructuredBlockForest( targetNumRootBlocks );
}
} // namespace mesh } // namespace mesh
} // namespace walberla } // namespace walberla
\ No newline at end of file
...@@ -102,7 +102,7 @@ void test( const std::string & meshFile, const uint_t numProcesses, const uint_t ...@@ -102,7 +102,7 @@ void test( const std::string & meshFile, const uint_t numProcesses, const uint_t
WALBERLA_LOG_INFO_ON_ROOT( "Creating SBF with StaticLevelwiseCurveBalanceWeighted Partitioner" ); WALBERLA_LOG_INFO_ON_ROOT( "Creating SBF with StaticLevelwiseCurveBalanceWeighted Partitioner" );
bfc.setTargetProcessAssignmentFunction( blockforest::StaticLevelwiseCurveBalanceWeighted() ); bfc.setTargetProcessAssignmentFunction( blockforest::StaticLevelwiseCurveBalanceWeighted() );
auto sbf_default = bfc.createSetupBlockForest( Vector3<uint_t>(64,64,64), numProcesses ); auto sbf_default = bfc.createSetupBlockForest( Vector3<uint_t>(64,64,64), numProcesses );
sbf_default->writeVTKOutput("sbf_default"); //sbf_default->writeVTKOutput("sbf_default");
WALBERLA_LOG_INFO_ON_ROOT( sbf_default->toString() ); WALBERLA_LOG_INFO_ON_ROOT( sbf_default->toString() );
return; return;
...@@ -137,6 +137,37 @@ void test( const std::string & meshFile, const uint_t numProcesses, const uint_t ...@@ -137,6 +137,37 @@ void test( const std::string & meshFile, const uint_t numProcesses, const uint_t
#endif #endif
} }
template< typename MeshType >
void testHelperFunctions( const std::string & meshFile, const uint_t numTotalBlocks )
{
auto mesh = make_shared<MeshType>();
mesh::readAndBroadcast( meshFile, *mesh);
auto triDist = make_shared< mesh::TriangleDistance<MeshType> >( mesh );
auto distanceOctree = make_shared< DistanceOctree< MeshType > >( triDist );
const real_t meshVolume = real_c( computeVolume( *mesh ) );
const real_t blockVolume = meshVolume / real_c( numTotalBlocks );
static const real_t cellsPersBlock = real_t(1000);
const real_t cellVolume = blockVolume / cellsPersBlock;
const real_t dx = std::pow( cellVolume, real_t(1) / real_t(3) );
WALBERLA_LOG_INFO_ON_ROOT( "Creating SBF with createStructuredBlockStorageInsideMesh with block size" );
auto sbf0 = mesh::createStructuredBlockStorageInsideMesh( distanceOctree, dx, numTotalBlocks );
WALBERLA_LOG_INFO_ON_ROOT( "Creating SBF with createStructuredBlockStorageInsideMesh with block size" );
Vector3<uint_t> blockSize( sbf0->getNumberOfXCells(), sbf0->getNumberOfYCells(), sbf0->getNumberOfZCells() );
auto sbf1 = mesh::createStructuredBlockStorageInsideMesh( distanceOctree, dx, blockSize );
auto exteriorAabb = computeAABB( *mesh ).getScaled( real_t(2) );
WALBERLA_LOG_INFO_ON_ROOT( "Creating SBF with createStructuredBlockStorageInsideMesh with block size" );
auto sbf2 = mesh::createStructuredBlockStorageOutsideMesh( exteriorAabb, distanceOctree, dx, numTotalBlocks );
WALBERLA_LOG_INFO_ON_ROOT( "Creating SBF with createStructuredBlockStorageInsideMesh with block size" );
blockSize = Vector3<uint_t>( sbf2->getNumberOfXCells(), sbf2->getNumberOfYCells(), sbf2->getNumberOfZCells() );
auto sbf3 = mesh::createStructuredBlockStorageOutsideMesh( exteriorAabb, distanceOctree, dx, blockSize );
}
int main( int argc, char * argv[] ) int main( int argc, char * argv[] )
{ {
debug::enterTestMode(); debug::enterTestMode();
...@@ -155,6 +186,8 @@ int main( int argc, char * argv[] ) ...@@ -155,6 +186,8 @@ int main( int argc, char * argv[] )
//test< mesh::FloatTriangleMesh >( meshFile, numProcesses, numTotalBlocks ); //test< mesh::FloatTriangleMesh >( meshFile, numProcesses, numTotalBlocks );
//test< mesh::PythonTriangleMesh >( meshFile, numProcesses, numTotalBlocks ); //test< mesh::PythonTriangleMesh >( meshFile, numProcesses, numTotalBlocks );
testHelperFunctions< mesh::TriangleMesh >( meshFile, numTotalBlocks );
return EXIT_SUCCESS; return EXIT_SUCCESS;
} }
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment