/* Lattice Boltzmann sample, written in C++, using the OpenLB
* library
*
* Copyright (C) 2007, 2012 Jonas Latt, Mathias J. Krause
* Vojtech Cvrcek, Peter Weisbrod
* E-mail contact: info@openlb.net
* The most recent release of OpenLB can be downloaded at
*
*
* This program 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 2
* of the License, or (at your option) any later version.
*
* This program 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 this program; if not, write to the Free
* Software Foundation, Inc., 51 Franklin Street, Fifth Floor,
* Boston, MA 02110-1301, USA.
*/
#include "olb2D.h"
#include "olb2D.hh"
#include
using namespace olb;
using namespace olb::descriptors;
typedef double T;
#define DESCRIPTOR D2Q9Descriptor
const T lx = 4.0; // length of the channel
const T ly = 1.0; // height of the channel
const int N = 30; // resolution of the model
const T Re = 100.; // Reynolds number
const T maxPhysT = 60.; // max. simulation time in s, SI unit
const T physInterval = 0.25; // interval for the convergence check in s
const T residuum = 1e-5; // residuum for the convergence check
void prepareGeometry(UnitConverter const& converter,
SuperGeometry2D& superGeometry)
{
OstreamManager clout(std::cout,"prepareGeometry");
clout << "Prepare Geometry ..." << std::endl;
superGeometry.rename(0,1);
const T physSpacing = converter.getPhysDeltaX();
// Set material number for bounce back boundaries
{
const Vector wallExtend {lx+physSpacing, physSpacing/2};
const Vector wallOrigin {-physSpacing/4, -physSpacing/4};
IndicatorCuboid2D lowerWall(wallExtend, wallOrigin);
superGeometry.rename(1,2,lowerWall);
IndicatorCuboid2D upperWall(wallExtend, wallOrigin + Vector {0,ly});
superGeometry.rename(1,2,upperWall);
}
// Set material number for inflow and outflow
{
const Vector extend { physSpacing/2, ly};
const Vector origin {-physSpacing/4, -physSpacing/4};
IndicatorCuboid2D inflow(extend, origin);
superGeometry.rename(1,3,inflow);
IndicatorCuboid2D outflow(extend, origin + Vector {lx,0});
superGeometry.rename(1,4,outflow);
}
// Set material number for vertically centered obstacle
{
const Vector extend {0.1, 0.1};
const Vector origin {1.25, (ly-extend[1])/2};
IndicatorCuboid2D obstacle(extend, origin);
superGeometry.rename(1,2,obstacle);
}
superGeometry.clean();
superGeometry.innerClean();
superGeometry.checkForErrors();
superGeometry.print();
clout << "Prepare Geometry ... OK" << std::endl;
}
void prepareCoarseLattice(UnitConverter const& converter,
SuperLattice2D& sLattice,
Dynamics& bulkDynamics,
sOnLatticeBoundaryCondition2D& sBoundaryCondition,
SuperGeometry2D& superGeometry)
{
OstreamManager clout(std::cout,"prepareLattice");
clout << "Prepare coarse lattice ..." << std::endl;
const T omega = converter.getLatticeRelaxationFrequency();
sLattice.defineDynamics(superGeometry, 0, &instances::getNoDynamics());
sLattice.defineDynamics(superGeometry, 1, &bulkDynamics); // bulk
sLattice.defineDynamics(superGeometry, 2, &instances::getBounceBack());
sLattice.defineDynamics(superGeometry, 3, &bulkDynamics); // inflow
sLattice.defineDynamics(superGeometry, 4, &bulkDynamics); // outflow
sBoundaryCondition.addVelocityBoundary(superGeometry, 3, omega);
sBoundaryCondition.addPressureBoundary(superGeometry, 4, omega);
const T Lx = converter.getLatticeLength(lx);
const T Ly = converter.getLatticeLength(ly);
const T p0 = 8.*converter.getLatticeViscosity()*converter.getCharLatticeVelocity()*Lx/(Ly*Ly);
AnalyticalLinear2D rho(-p0/lx*DESCRIPTOR::invCs2, 0, p0*DESCRIPTOR::invCs2+1);
const T maxVelocity = converter.getCharLatticeVelocity();
const T radius = ly/2;
std::vector axisPoint{0, ly/2};
std::vector axisDirection{1, 0};
Poiseuille2D u(axisPoint, axisDirection, maxVelocity, radius);
sLattice.defineRhoU(superGeometry, 1, rho, u);
sLattice.iniEquilibrium(superGeometry, 1, rho, u);
sLattice.defineRhoU(superGeometry, 2, rho, u);
sLattice.iniEquilibrium(superGeometry, 2, rho, u);
sLattice.defineRhoU(superGeometry, 3, rho, u);
sLattice.iniEquilibrium(superGeometry, 3, rho, u);
sLattice.defineRhoU(superGeometry, 4, rho, u);
sLattice.iniEquilibrium(superGeometry, 4, rho, u);
sLattice.initialize();
clout << "Prepare coarse lattice ... OK" << std::endl;
}
void prepareFineLattice(UnitConverter const& converter,
SuperLattice2D& sLattice,
Dynamics& bulkDynamics,
sOnLatticeBoundaryCondition2D& sBoundaryCondition,
SuperGeometry2D& superGeometry)
{
OstreamManager clout(std::cout,"prepareLattice");
clout << "Prepare fine lattice ..." << std::endl;
const T omega = converter.getLatticeRelaxationFrequency();
sLattice.defineDynamics(superGeometry, 0, &instances::getNoDynamics());
sLattice.defineDynamics(superGeometry, 1, &bulkDynamics); // bulk
sLattice.defineDynamics(superGeometry, 2, &instances::getBounceBack());
const T Lx = converter.getLatticeLength(lx);
const T Ly = converter.getLatticeLength(ly);
const T p0 = 8.*converter.getLatticeViscosity()*converter.getCharLatticeVelocity()*Lx/(Ly*Ly);
AnalyticalLinear2D rho(-p0/lx*DESCRIPTOR::invCs2, 0, p0*DESCRIPTOR::invCs2+1);
const T maxVelocity = converter.getCharLatticeVelocity();
const T radius = ly/2;
std::vector axisPoint{0, ly/2};
std::vector axisDirection{1, 0};
Poiseuille2D u(axisPoint, axisDirection, maxVelocity, radius);
sLattice.defineRhoU(superGeometry, 1, rho, u);
sLattice.iniEquilibrium(superGeometry, 1, rho, u);
sLattice.defineRhoU(superGeometry, 2, rho, u);
sLattice.iniEquilibrium(superGeometry, 2, rho, u);
sLattice.initialize();
clout << "Prepare fine lattice ... OK" << std::endl;
}
void getResults(const std::string& prefix,
SuperLattice2D& sLattice,
UnitConverter const& converter, int iT,
SuperGeometry2D& superGeometry, Timer& timer, bool hasConverged)
{
OstreamManager clout(std::cout,"getResults");
SuperVTMwriter2D vtmWriter(prefix + "poiseuille2d");
SuperLatticePhysVelocity2D velocity(sLattice, converter);
SuperLatticePhysPressure2D pressure(sLattice, converter);
SuperLatticeGeometry2D geometry( sLattice, superGeometry );
vtmWriter.addFunctor(geometry);
vtmWriter.addFunctor(velocity);
vtmWriter.addFunctor(pressure);
const int statIter = converter.getLatticeTime(maxPhysT/10.);
if (iT==0) {
vtmWriter.createMasterFile();
}
if (iT%100==0) {
vtmWriter.write(iT);
}
if (iT%statIter==0 || hasConverged) {
timer.update(iT);
timer.printStep();
sLattice.getStatistics().print(iT,converter.getPhysTime(iT));
}
}
int main(int argc, char* argv[])
{
olbInit(&argc, &argv);
singleton::directories().setOutputDir("./tmp/");
OstreamManager clout(std::cout,"main");
const Vector coarseOrigin {0.0, 0.0};
const Vector coarseExtend {lx, ly};
IndicatorCuboid2D coarseCuboid(coarseExtend, coarseOrigin);
auto coarseGrid = Grid2D::make(coarseCuboid, N, 0.8, Re);
prepareGeometry(coarseGrid->getConverter(), coarseGrid->getSuperGeometry());
const Vector wantedFineExtend {2.0, 0.75};
const Vector fineOrigin = coarseGrid->alignToGrid({0.8, (ly-wantedFineExtend[1])/2});
const Vector fineExtend = coarseGrid->alignToGrid(fineOrigin + wantedFineExtend) - fineOrigin;
auto fineGrid = &coarseGrid->refine(fineOrigin, fineExtend);
prepareGeometry(fineGrid->getConverter(), fineGrid->getSuperGeometry());
const T coarseDeltaX = coarseGrid->getConverter().getPhysDeltaX();
IndicatorCuboid2D refinedOverlap(fineExtend - 4*coarseDeltaX, fineOrigin + 2*coarseDeltaX);
coarseGrid->getSuperGeometry().rename(1,0,refinedOverlap);
coarseGrid->getSuperGeometry().rename(2,0,refinedOverlap);
Dynamics* coarseBulkDynamics;
coarseBulkDynamics = new BGKdynamics(
coarseGrid->getConverter().getLatticeRelaxationFrequency(),
instances::getBulkMomenta());
sOnLatticeBoundaryCondition2D coarseSBoundaryCondition(coarseGrid->getSuperLattice());
createLocalBoundaryCondition2D(coarseSBoundaryCondition);
prepareCoarseLattice(
coarseGrid->getConverter(),
coarseGrid->getSuperLattice(),
*coarseBulkDynamics,
coarseSBoundaryCondition,
coarseGrid->getSuperGeometry());
Dynamics* fineBulkDynamics;
fineBulkDynamics = new BGKdynamics(
fineGrid->getConverter().getLatticeRelaxationFrequency(),
instances::getBulkMomenta());
sOnLatticeBoundaryCondition2D fineSBoundaryCondition(fineGrid->getSuperLattice());
createLocalBoundaryCondition2D(fineSBoundaryCondition);
prepareFineLattice(
fineGrid->getConverter(),
fineGrid->getSuperLattice(),
*fineBulkDynamics,
fineSBoundaryCondition,
fineGrid->getSuperGeometry());
clout << "starting simulation..." << endl;
Timer timer(
coarseGrid->getConverter().getLatticeTime(maxPhysT),
coarseGrid->getSuperGeometry().getStatistics().getNvoxel());
util::ValueTracer converge(
fineGrid->getConverter().getLatticeTime(physInterval),
residuum);
timer.start();
for (int iT = 0; iT < coarseGrid->getConverter().getLatticeTime(maxPhysT); ++iT) {
if (converge.hasConverged()) {
clout << "Simulation converged." << endl;
break;
}
coarseGrid->collideAndStream();
getResults(
"coarse_",
coarseGrid->getSuperLattice(),
coarseGrid->getConverter(),
iT,
coarseGrid->getSuperGeometry(),
timer,
converge.hasConverged());
getResults(
"fine_",
fineGrid->getSuperLattice(),
fineGrid->getConverter(),
iT,
fineGrid->getSuperGeometry(),
timer,
converge.hasConverged());
converge.takeValue(fineGrid->getSuperLattice().getStatistics().getAverageEnergy(), true);
}
timer.stop();
timer.printSummary();
}