diff --git a/apps/tutorials/CMakeLists.txt b/apps/tutorials/CMakeLists.txt index 1e4a57916efb9b93d0fe3f47fdd1395bc83282ad..91dae24cd697064fe40b1490b45f7a863f4fffab 100644 --- a/apps/tutorials/CMakeLists.txt +++ b/apps/tutorials/CMakeLists.txt @@ -1,4 +1,4 @@ add_subdirectory(basics) add_subdirectory(lbm) add_subdirectory(pde) - \ No newline at end of file +add_subdirectory(pe) diff --git a/apps/tutorials/pe/01_ConfinedGas.cpp b/apps/tutorials/pe/01_ConfinedGas.cpp new file mode 100644 index 0000000000000000000000000000000000000000..4f0190197ff14bc26b22daf7ad2afb35dc530889 --- /dev/null +++ b/apps/tutorials/pe/01_ConfinedGas.cpp @@ -0,0 +1,161 @@ +//====================================================================================================================== +// +// This file is part of waLBerla. waLBerla is free software: you can +// redistribute it and/or modify it under the terms of the GNU General Public +// License as published by the Free Software Foundation, either version 3 of +// the License, or (at your option) any later version. +// +// waLBerla is distributed in the hope that it will be useful, but WITHOUT +// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or +// FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License +// for more details. +// +// You should have received a copy of the GNU General Public License along +// with waLBerla (see COPYING.txt). If not, see <http://www.gnu.org/licenses/>. +// +//! \file 01_ConfinedGas.cpp +//! \author Sebastian Eibl <sebastian.eibl@fau.de> +// +//====================================================================================================================== + +//! [Includes] +#include <pe/basic.h> + +#include <core/Environment.h> +#include <core/grid_generator/HCPIterator.h> +#include <core/grid_generator/SCIterator.h> +#include <core/logging/Logging.h> +//! [Includes] + +using namespace walberla; +using namespace walberla::pe; + +//! [BodyTypeTuple] +typedef boost::tuple<Sphere, Plane> BodyTypeTuple ; +//! [BodyTypeTuple] + +int main( int argc, char ** argv ) +{ + //! [Parameters] + Environment env(argc, argv); + WALBERLA_UNUSED(env); + + math::seedRandomGenerator( static_cast<unsigned int>(1337 * mpi::MPIManager::instance()->worldRank()) ); + + real_t spacing = real_c(1.0); + real_t radius = real_c(0.4); + real_t vMax = real_c(1.0); + int simulationSteps = 200; + real_t dt = real_c(0.01); + //! [Parameters] + + WALBERLA_LOG_INFO("*** GLOBALBODYSTORAGE ***"); + //! [GlobalBodyStorage] + shared_ptr<BodyStorage> globalBodyStorage = make_shared<BodyStorage>(); + //! [GlobalBodyStorage] + + WALBERLA_LOG_INFO("*** BLOCKFOREST ***"); + // create forest + //! [BlockForest] + shared_ptr< BlockForest > forest = createBlockForest( AABB(0,0,0,20,20,20), // simulation domain + Vector3<uint_t>(1,1,1), // blocks in each direction + Vector3<bool>(false, false, false) // periodicity + ); + //! [BlockForest] + if (!forest) + { + WALBERLA_LOG_INFO( "No BlockForest created ... exiting!"); + return EXIT_SUCCESS; + } + + WALBERLA_LOG_INFO("*** STORAGEDATAHANDLING ***"); + // add block data + //! [StorageDataHandling] + auto storageID = forest->addBlockData(createStorageDataHandling<BodyTypeTuple>(), "Storage"); + //! [StorageDataHandling] + //! [AdditionalBlockData] + auto ccdID = forest->addBlockData(ccd::createHashGridsDataHandling( globalBodyStorage, storageID ), "CCD"); + auto fcdID = forest->addBlockData(fcd::createSimpleFCDDataHandling<BodyTypeTuple>(), "FCD"); + //! [AdditionalBlockData] + + WALBERLA_LOG_INFO("*** INTEGRATOR ***"); + //! [Integrator] + cr::HCSITS cr(globalBodyStorage, forest, storageID, ccdID, fcdID); + cr.setMaxIterations( 10 ); + cr.setRelaxationModel( cr::HardContactSemiImplicitTimesteppingSolvers::ApproximateInelasticCoulombContactByDecoupling ); + cr.setRelaxationParameter( real_t(0.7) ); + cr.setGlobalLinearAcceleration( Vec3(0,0,0) ); + //! [Integrator] + + WALBERLA_LOG_INFO("*** BodyTypeTuple ***"); + // initialize body type ids + //! [SetBodyTypeIDs] + SetBodyTypeIDs<BodyTypeTuple>::execute(); + //! [SetBodyTypeIDs] + + WALBERLA_LOG_INFO("*** SETUP - START ***"); + //! [Material] + const real_t static_cof ( real_c(0.1) / 2 ); // Coefficient of static friction. Note: pe doubles the input coefficient of friction for material-material contacts. + const real_t dynamic_cof ( static_cof ); // Coefficient of dynamic friction. Similar to static friction for low speed friction. + MaterialID material = createMaterial( "granular", real_t( 1.0 ), 0, static_cof, dynamic_cof, real_t( 0.5 ), 1, 1, 0, 0 ); + //! [Material] + + auto simulationDomain = forest->getDomain(); + auto generationDomain = simulationDomain; // simulationDomain.getExtended(-real_c(0.5) * spacing); + //! [Planes] + createPlane(*globalBodyStorage, 0, Vec3(1,0,0), simulationDomain.minCorner(), material ); + createPlane(*globalBodyStorage, 0, Vec3(-1,0,0), simulationDomain.maxCorner(), material ); + createPlane(*globalBodyStorage, 0, Vec3(0,1,0), simulationDomain.minCorner(), material ); + createPlane(*globalBodyStorage, 0, Vec3(0,-1,0), simulationDomain.maxCorner(), material ); + createPlane(*globalBodyStorage, 0, Vec3(0,0,1), simulationDomain.minCorner(), material ); + createPlane(*globalBodyStorage, 0, Vec3(0,0,-1), simulationDomain.maxCorner(), material ); + //! [Planes] + + //! [Gas] + uint_t numParticles = uint_c(0); + for (auto blkIt = forest->begin(); blkIt != forest->end(); ++blkIt) + { + IBlock & currentBlock = *blkIt; + for (auto it = grid_generator::SCIterator(currentBlock.getAABB().getIntersection(generationDomain), Vector3<real_t>(spacing, spacing, spacing) * real_c(0.5), spacing); it != grid_generator::SCIterator(); ++it) + { + SphereID sp = createSphere( *globalBodyStorage, *forest, storageID, 0, *it, radius, material); + Vec3 rndVel(math::realRandom<real_t>(-vMax, vMax), math::realRandom<real_t>(-vMax, vMax), math::realRandom<real_t>(-vMax, vMax)); + if (sp != NULL) sp->setLinearVel(rndVel); + if (sp != NULL) ++numParticles; + } + } + WALBERLA_LOG_INFO("#particles created: " << numParticles); + //! [Gas] + + WALBERLA_LOG_INFO("*** SETUP - END ***"); + + WALBERLA_LOG_INFO("*** SIMULATION - START ***"); + //! [GameLoop] + for (int i=0; i < simulationSteps; ++i) + { + if( i % 10 == 0 ) + { + WALBERLA_LOG_DEVEL( "Timestep " << i << " / " << simulationSteps ); + } + + cr.timestep( real_c(dt) ); + } + //! [GameLoop] + WALBERLA_LOG_INFO("*** SIMULATION - END ***"); + + WALBERLA_LOG_INFO("*** GETTING STATISTICAL INFORMATION ***"); + //! [PostProcessing] + Vec3 meanVelocity(0,0,0); + for (auto blockIt = forest->begin(); blockIt != forest->end(); ++blockIt) + { + for (auto bodyIt = LocalBodyIterator::begin(*blockIt, storageID); bodyIt != LocalBodyIterator::end(); ++bodyIt) + { + meanVelocity += bodyIt->getLinearVel(); + } + } + meanVelocity /= numParticles; + WALBERLA_LOG_INFO( "mean velocity: " << meanVelocity ); + //! [PostProcessing] + + return EXIT_SUCCESS; +} diff --git a/apps/tutorials/pe/01_ConfinedGas.dox b/apps/tutorials/pe/01_ConfinedGas.dox new file mode 100644 index 0000000000000000000000000000000000000000..28fc0a834f5b6a70d393e94a0c989be6ab04ee1d --- /dev/null +++ b/apps/tutorials/pe/01_ConfinedGas.dox @@ -0,0 +1,78 @@ +namespace walberla { +namespace pe { + +/** +\page tutorial_pe_01 Tutorial - Confined Gas: Setting up a simple pe simulation + +This tutorial will help you to setup up a simple pe simulation. The scenario for this tutorial a is particle +gas (spheres) confined by solid walls. + +For this tutorial only a few includes are needed, namely: +\snippet 01_ConfinedGas.cpp Includes + +The first important step is to define a BodyTuple. This tuple stores all body types which will be used +during the simulation. This tuple helps to generate all the functions which need to be specialized for +a particular body type. +\attention If you miss a body type the behaviour of the simulation is undefined. Most probably there will be an error. + +\snippet 01_ConfinedGas.cpp BodyTypeTuple + +Next the waLBerla environment is initalized, the random number generator is seeded and some simulation parameters are set. +\snippet 01_ConfinedGas.cpp Parameters + +The BlockForest is the main datastructure in the waLBerla framework. It is responsible for the domain decomposition and +holds all the blocks with their data. For more information about the general design of the waLBerla framework please refer +to \ref tutorial_basics_01 and the documentation of domain_decomposition::BlockStorage. For this tutorial the number of blocks +in each direction is fixed to 1 and has to stay like that. In tutorial 2 we will talk about using more than one block and +parallelism. +\snippet 01_ConfinedGas.cpp BlockForest + +There are two types of storages to store particle information. One is the global body storage which is responsible for very +large particles. These particles are stored on all processes. +\snippet 01_ConfinedGas.cpp GlobalBodyStorage +The other one is a block local storage for all bodies which are located within the subdomain of a block. You can find more +information about this concept of data on blocks in \ref tutorial_basics_01. +\snippet 01_ConfinedGas.cpp StorageDataHandling + +In addition to the block local storage also the coarse as well as the fine collision detection needs storage on each block. +Therefore the corresponding data handling has to be registered. +\snippet 01_ConfinedGas.cpp AdditionalBlockData + +Only one final component is missing for a successfull simulation - the time integrator. Currently there exists two +integrators cr::DEM for soft contacts and cr::HCSITS for hard contacts. These have to be setup as local objects. +\snippet 01_ConfinedGas.cpp Integrator + +We now have all our tools together and can start creating the simulation. First we have to specify the material parameters +we want to use. +\snippet 01_ConfinedGas.cpp Material + +And we have to initialize all body type ids. +\snippet 01_ConfinedGas.cpp SetBodyTypeIDs + +Then we can set up the confining walls. +\snippet 01_ConfinedGas.cpp Planes + +The gas particles are generated with the help of grid generators. These generators are iterators which facilitate the +construction of particles on a regular grid. Currently grid_generator::SCIterator (for a simple cubic lattice) as well +as grid_generator::HCPIterator (for a hexagonal close packing lattice) are supported. The spheres are create with createSphere +which returns a SphereID of the created sphere. This SphereID acts like a pointer to Sphere and can be used like that. +If for various reasons the sphere was not created the return value is NULL. +\attention Before accessing the underlying sphere you should always check for NULL! + +\snippet 01_ConfinedGas.cpp Gas + +Since the setup is finished now we can run the simulation loop. The simulation loop is as simple as: +\snippet 01_ConfinedGas.cpp GameLoop + +After the simulation is finished we can collect the results. In this case we only calculate the mean velocity of all particles. +The particles can be easily accessed via the LocalBodyIterator. This iterator allows us to iterate through all particles and +access their properties. +\snippet 01_ConfinedGas.cpp PostProcessing + +Congratulation! You successfully created your first rigid body simulation. +In the next tutorial we will look at possible extensions which can make your live easier as well as how to run parallel simulations. + +*/ + +} +} diff --git a/apps/tutorials/pe/CMakeLists.txt b/apps/tutorials/pe/CMakeLists.txt new file mode 100644 index 0000000000000000000000000000000000000000..1a236f44933c947cada7ee5a84c4a8e2a698063b --- /dev/null +++ b/apps/tutorials/pe/CMakeLists.txt @@ -0,0 +1,5 @@ +waLBerla_link_files_to_builddir( *.cfg ) + +waLBerla_add_executable ( NAME 01_Tutorial_ConfinedGas + FILES 01_ConfinedGas.cpp + DEPENDS blockforest core pe ) diff --git a/doc/Mainpage.dox b/doc/Mainpage.dox index b31a3a182fe3934b774ff486b71e419753505dc6..e89de056aaae2f09261def5c76c827e6eea398ca 100644 --- a/doc/Mainpage.dox +++ b/doc/Mainpage.dox @@ -16,11 +16,13 @@ all the basic data strcutures and concepts of the framework. - \ref tutorial_basics_01 \n The first tutorial explains the basic data structures of waLBerla and how to setup a simple domain. - \ref tutorial_basics_02 \n - The second tutorial shows how to write algorithms that operate on block data. + The second tutorial shows how to write algorithms that operate on block data. - \ref tutorial_basics_03 \n - The third tutorial deals with communication, iterating field/grid/lattice data, and geometry setup. + The third tutorial deals with communication, iterating field/grid/lattice data, and geometry setup. +\subsection pe Rigid Body Dynamics (pe) +- \ref tutorial_pe_01 \n \subsection pdes Solving Partial Differential Equations @@ -33,11 +35,6 @@ all the basic data strcutures and concepts of the framework. - \ref tutorial_lbm01 \n A full LBM simulation is built. - -\section TechDetails Technical Details - -- \ref TechDetailsPe \n - Technical Details about the pe physics module. \section cite Please cite us