Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
//======================================================================================================================
//
// 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 03_AdvancedLBMCodegen.cpp
//! \author Frederik Hennig <frederik.hennig@fau.de>
//
//======================================================================================================================
#include "blockforest/all.h"
#include "core/all.h"
#include "domain_decomposition/all.h"
#include "field/all.h"
#include "geometry/all.h"
#include "lbm/boundary/factories/DefaultBoundaryHandling.h"
#include "lbm/communication/PdfFieldPackInfo.h"
#include "lbm/field/AddToStorage.h"
#include "lbm/field/PdfField.h"
#include "lbm/field/initializer/all.h"
#include "lbm/vtk/VTKOutput.h"
#include "stencil/D2Q9.h"
#include "timeloop/all.h"
// Codegen Includes
#include "CumulantMRTNoSlip.h"
#include "CumulantMRTPackInfo.h"
#include "CumulantMRTSweep.h"
#include "DensityAndVelocityFieldSetter.h"
namespace walberla
{
///////////////////////
/// Typedef Aliases ///
///////////////////////
// Communication Pack Info
typedef pystencils::CumulantMRTPackInfo PackInfo_T;
// LB Method Stencil
typedef stencil::D2Q9 Stencil_T;
// PDF field type
typedef field::GhostLayerField< real_t, Stencil_T::Size > PdfField_T;
// Velocity Field Type
typedef field::GhostLayerField< real_t, Stencil_T::D > VectorField_T;
// Boundary Handling
typedef walberla::uint8_t flag_t;
typedef FlagField< flag_t > FlagField_T;
typedef lbm::CumulantMRTNoSlip NoSlip_T;
//////////////////////////////////////////
/// Shear Flow Velocity Initialization ///
//////////////////////////////////////////
bool initShearFlowVelocityField(const shared_ptr< StructuredBlockForest >& blocks, const BlockDataID& velocityFieldId,
const Config::BlockHandle& config)
{
math::RealRandom< real_t > rng(config.getParameter< std::mt19937::result_type >("noise_seed", 42));
real_t velocityMagnitude = config.getParameter< real_t >("VelocityMagnitude", real_c(0.08));
real_t noiseMagnitude = config.getParameter< real_t >("NoiseMagnitude", real_c(0.1) * velocityMagnitude);
real_t n_y = real_c(blocks->getNumberOfYCells());
for (auto blockIt = blocks->begin(); blockIt != blocks->end(); ++blockIt)
{
auto u = (*blockIt).getData< VectorField_T >(velocityFieldId);
for (auto cellIt = u->beginWithGhostLayerXYZ(); cellIt != u->end(); ++cellIt)
{
Cell globalCell(cellIt.cell());
blocks->transformBlockLocalToGlobalCell(globalCell, *blockIt);
real_t relative_y = real_c(globalCell.y()) / n_y;
u->get(cellIt.cell(), 0) = relative_y < 0.3 || relative_y > 0.7 ? velocityMagnitude : -velocityMagnitude;
u->get(cellIt.cell(), 1) = noiseMagnitude * rng();
}
}
}
/////////////////////
/// Main Function ///
/////////////////////
int main(int argc, char** argv)
{
walberla::Environment walberlaEnv(argc, argv);
if (!walberlaEnv.config()) { WALBERLA_ABORT("No configuration file specified!"); }
///////////////////////////////////////////////////////
/// Block Storage Creation and Simulation Parameter ///
///////////////////////////////////////////////////////
auto blocks = blockforest::createUniformBlockGridFromConfig(walberlaEnv.config());
// read parameters
auto parameters = walberlaEnv.config()->getOneBlock("Parameters");
const uint_t timesteps = parameters.getParameter< uint_t >("timesteps", uint_c(10));
const real_t omega = parameters.getParameter< real_t >("omega", real_c(1.8));
const double remainingTimeLoggerFrequency =
parameters.getParameter< double >("remainingTimeLoggerFrequency", 3.0); // in seconds
////////////////////////////
/// Velocity Field Setup ///
////////////////////////////
BlockDataID velocityFieldId = field::addToStorage< VectorField_T >(blocks, "velocity", real_c(0.0), field::fzyx);
///////////////////////
/// PDF Field Setup ///
///////////////////////
BlockDataID pdfFieldId = field::addToStorage< PdfField_T >(blocks, "pdf field", real_c(0.0), field::fzyx);
BlockDataID flagFieldId = field::addFlagFieldToStorage< FlagField_T >(blocks, "flag field");
// First, set up the initial shear flow in the velocity field...
auto shearFlowSetup = walberlaEnv.config()->getOneBlock("ShearFlowSetup");
initShearFlowVelocityField(blocks, velocityFieldId, shearFlowSetup);
real_t rho = shearFlowSetup.getParameter("rho", real_c(1.8));
// ... and then use the generated PDF setter to initialize the PDF field.
pystencils::DensityAndVelocityFieldSetter pdfSetter(pdfFieldId, velocityFieldId, rho);
for(auto blockIt = blocks->begin(); blockIt != blocks->end(); ++blockIt){
pdfSetter(&(*blockIt));
}
/////////////////////////
/// Boundary Handling ///
/////////////////////////
const FlagUID fluidFlagUID("Fluid");
auto boundariesConfig = walberlaEnv.config()->getOneBlock("Boundaries");
// BlockDataID boundaryHandlingId =
// BHFactory::addBoundaryHandlingToStorage(blocks, "boundary handling", flagFieldId, pdfFieldId, fluidFlagUID,
// Vector3< real_t >(), Vector3< real_t >(), real_c(0.0), real_c(0.0));
// geometry::initBoundaryHandling< BHFactory::BoundaryHandling >(*blocks, boundaryHandlingId, boundariesConfig);
// geometry::setNonBoundaryCellsToDomain< BHFactory::BoundaryHandling >(*blocks, boundaryHandlingId);
/////////////////
/// Time Loop ///
/////////////////
SweepTimeloop timeloop(blocks->getBlockStorage(), timesteps);
// Communication
blockforest::communication::UniformBufferedScheme< Stencil_T > communication(blocks);
communication.addPackInfo(make_shared< PackInfo_T >(pdfFieldId));
// Timeloop
timeloop.add() << BeforeFunction(communication, "communication");
// Stability Checker
timeloop.addFuncAfterTimeStep(makeSharedFunctor(field::makeStabilityChecker< PdfField_T, FlagField_T >(
walberlaEnv.config(), blocks, pdfFieldId, flagFieldId, fluidFlagUID)),
"LBM stability check");
// Time logger
timeloop.addFuncAfterTimeStep(timing::RemainingTimeLogger(timeloop.getNrOfTimeSteps(), remainingTimeLoggerFrequency),
"remaining time logger");
// TODO: Replace
// LBM VTK Output
// lbm::VTKOutput< LatticeModel_T, FlagField_T >::addToTimeloop(timeloop, blocks, walberlaEnv.config(), pdfFieldId,
// flagFieldId, fluidFlagUID);
timeloop.run();
return EXIT_SUCCESS;
}
} // namespace walberla
int main(int argc, char** argv) { return walberla::main(argc, argv); }