/*
* Lattice Boltzmann grid refinement sample, written in C++,
* using the OpenLB library
*
* Copyright (C) 2019 Adrian Kummerländer
* 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"
#ifndef OLB_PRECOMPILED
#include "olb2D.hh"
#endif
#include
using namespace olb;
typedef double T;
#define DESCRIPTOR descriptors::D2Q9Descriptor
const T lx = 4.0; // length of the channel
const T ly = 1.0; // height of the channel
const int N = 10; // resolution of the model
const T Re = 10.; // Reynolds number
const T uMax = 0.01; // Max lattice speed
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(Grid2D& grid)
{
OstreamManager clout(std::cout,"prepareGeometry");
clout << "Prepare Geometry ..." << std::endl;
auto& converter = grid.getConverter();
auto& sGeometry = grid.getSuperGeometry();
sGeometry.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);
sGeometry.rename(1,2,lowerWall);
IndicatorCuboid2D upperWall(wallExtend, wallOrigin + Vector {0,ly});
sGeometry.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);
sGeometry.rename(1,3,inflow);
IndicatorCuboid2D outflow(extend, origin + Vector {lx,0});
sGeometry.rename(1,4,outflow);
}
sGeometry.clean();
sGeometry.innerClean();
sGeometry.checkForErrors();
sGeometry.print();
clout << "Prepare Geometry ... OK" << std::endl;
}
void prepareLattice(Grid2D& grid)
{
OstreamManager clout(std::cout,"prepareLattice");
clout << "Prepare lattice ..." << std::endl;
auto& converter = grid.getConverter();
auto& sGeometry = grid.getSuperGeometry();
auto& sLattice = grid.getSuperLattice();
Dynamics& bulkDynamics = grid.addDynamics(
std::unique_ptr>(
new BGKdynamics(
grid.getConverter().getLatticeRelaxationFrequency(),
instances::getBulkMomenta())));
sOnLatticeBoundaryCondition2D& sBoundaryCondition = grid.getOnLatticeBoundaryCondition();
createInterpBoundaryCondition2D(sBoundaryCondition);
const T omega = converter.getLatticeRelaxationFrequency();
sLattice.defineDynamics(sGeometry, 0, &instances::getNoDynamics());
sLattice.defineDynamics(sGeometry, 1, &bulkDynamics); // bulk
sLattice.defineDynamics(sGeometry, 2, &bulkDynamics); // walls
sLattice.defineDynamics(sGeometry, 3, &bulkDynamics); // inflow
sLattice.defineDynamics(sGeometry, 4, &bulkDynamics); // outflow
sBoundaryCondition.addVelocityBoundary(sGeometry, 2, omega); // 0-velocity walls
sBoundaryCondition.addVelocityBoundary(sGeometry, 3, omega); // velocity inflow
sBoundaryCondition.addPressureBoundary(sGeometry, 4, omega); // pressure outflow
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);
auto materials = sGeometry.getMaterialIndicator({1,2,3,4});
sLattice.defineRhoU(materials, rho, u);
sLattice.iniEquilibrium(materials, rho, u);
sLattice.initialize();
clout << "Prepare lattice ... OK" << std::endl;
}
void getResults(const std::string& prefix,
Grid2D& grid,
int iT, Timer& timer, bool hasConverged)
{
OstreamManager clout(std::cout,"getResults");
auto& converter = grid.getConverter();
auto& sLattice = grid.getSuperLattice();
auto& sGeometry = grid.getSuperGeometry();
SuperVTMwriter2D vtmWriter(prefix + "poiseuille2d");
SuperLatticePhysVelocity2D velocity(sLattice, converter);
SuperLatticePhysPressure2D pressure(sLattice, converter);
SuperLatticeGeometry2D geometry(sLattice, sGeometry);
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));
}
}
void getError(const std::string& prefix,
Grid2D& grid)
{
OstreamManager clout(std::cout, prefix);
auto& converter = grid.getConverter();
auto& sLattice = grid.getSuperLattice();
auto& sGeometry = grid.getSuperGeometry();
const T maxVelocity = converter.getCharPhysVelocity();
const T radius = ly/2;
std::vector axisPoint{lx/2, ly/2};
std::vector axisDirection{1, 0};
Poiseuille2D uSol(axisPoint, axisDirection, maxVelocity, radius);
const T Lx = converter.getLatticeLength(lx);
const T Ly = converter.getLatticeLength(ly);
T p0 = 8.*converter.getLatticeViscosity()*converter.getCharLatticeVelocity()*Lx/T(Ly*Ly);
AnalyticalLinear2D pSol(-converter.getPhysPressure(p0)/lx, 0, converter.getPhysPressure(p0));
SuperLatticePhysVelocity2D u(sLattice, converter);
SuperLatticePhysPressure2D p(sLattice, converter);
auto fluid = sGeometry.getMaterialIndicator(1);
int tmp[]= { };
T result[2]= { };
SuperRelativeErrorL2Norm2D velocityError(u, uSol, fluid);
velocityError(result, tmp);
clout << "velocity-L2-error(rel)=" << result[0] << std::endl;
SuperRelativeErrorL2Norm2D pressureError(p, pSol, fluid);
pressureError(result, tmp);
clout << "pressure-L2-error(rel)=" << result[0] << std::endl;
}
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/2, ly};
IndicatorCuboid2D coarseCuboid(coarseExtend, coarseOrigin);
Grid2D coarseGrid(coarseCuboid, LatticeVelocity(uMax), N, Re);
prepareGeometry(coarseGrid);
const T coarseDeltaX = coarseGrid.getConverter().getPhysDeltaX();
const Vector fineExtend {lx/2+coarseDeltaX, ly};
const Vector fineOrigin {lx/2-coarseDeltaX, 0.0};
auto& fineGrid = coarseGrid.refine(fineOrigin, fineExtend, false);
{
const Vector origin = fineGrid.getOrigin();
const Vector extend = fineGrid.getExtend();
const Vector extendY = {0,extend[1]};
coarseGrid.addFineCoupling(fineGrid, origin, extendY);
coarseGrid.addCoarseCoupling(fineGrid, origin + Vector {coarseDeltaX,0}, extendY);
}
prepareGeometry(fineGrid);
prepareLattice(coarseGrid);
prepareLattice(fineGrid);
clout << "Total number of active cells: " << coarseGrid.getActiveVoxelN() << endl;
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,
iT,
timer,
converge.hasConverged());
getResults(
"fine_",
fineGrid,
iT,
timer,
converge.hasConverged());
converge.takeValue(fineGrid.getSuperLattice().getStatistics().getAverageEnergy(), true);
}
getError("coarse", coarseGrid);
getError("fine", fineGrid);
timer.stop();
timer.printSummary();
}