Skip to content
Projects
Groups
Snippets
Help
Loading...
Help
Support
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
W
waLBerla
Project overview
Project overview
Details
Activity
Releases
Cycle Analytics
Repository
Repository
Files
Commits
Branches
Tags
Contributors
Graph
Compare
Charts
Issues
22
Issues
22
List
Boards
Labels
Milestones
Merge Requests
4
Merge Requests
4
CI / CD
CI / CD
Pipelines
Jobs
Schedules
Charts
Wiki
Wiki
Members
Members
Collapse sidebar
Close sidebar
Activity
Graph
Charts
Create a new issue
Jobs
Commits
Issue Boards
Open sidebar
waLBerla
waLBerla
Commits
eb1fd933
Commit
eb1fd933
authored
Aug 12, 2019
by
Martin Bauer
Browse files
Options
Browse Files
Download
Plain Diff
Merge branch 'fluctuating_mrt' into 'master'
Fluctuating MRT Closes
#80
See merge request
!220
parents
4f97ea03
40c2b271
Pipeline
#17166
passed with stages
in 425 minutes and 26 seconds
Changes
3
Pipelines
2
Hide whitespace changes
Inline
Side-by-side
Showing
3 changed files
with
145 additions
and
0 deletions
+145
-0
tests/lbm/CMakeLists.txt
tests/lbm/CMakeLists.txt
+4
-0
tests/lbm/codegen/FluctuatingMRT.cpp
tests/lbm/codegen/FluctuatingMRT.cpp
+119
-0
tests/lbm/codegen/FluctuatingMRT.py
tests/lbm/codegen/FluctuatingMRT.py
+22
-0
No files found.
tests/lbm/CMakeLists.txt
View file @
eb1fd933
...
...
@@ -77,3 +77,7 @@ waLBerla_python_file_generates(codegen/LbCodeGenerationExample.py
LbCodeGenerationExample_UBB.cpp
)
waLBerla_compile_test
(
FILES codegen/LbCodeGenerationExample.py
codegen/LbCodeGenerationExample.cpp
)
waLBerla_python_file_generates
(
codegen/FluctuatingMRT.py
FluctuatingMRT_LatticeModel.cpp
)
waLBerla_compile_test
(
FILES codegen/FluctuatingMRT.py
codegen/FluctuatingMRT.cpp
)
\ No newline at end of file
tests/lbm/codegen/FluctuatingMRT.cpp
0 → 100644
View file @
eb1fd933
//======================================================================================================================
//
// 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 SrtWithForceField.cpp
//! \author Michael Kuron <mkuron@icp.uni-stuttgart.de>
//
//======================================================================================================================
#include "blockforest/all.h"
#include "core/all.h"
#include "domain_decomposition/all.h"
#include "field/all.h"
#include "geometry/all.h"
#include "gui/all.h"
#include "timeloop/all.h"
#include "lbm/field/PdfField.h"
#include "lbm/field/AddToStorage.h"
#include "lbm/communication/PdfFieldPackInfo.h"
#include "lbm/gui/Connection.h"
#include "lbm/vtk/VTKOutput.h"
#include "FluctuatingMRT_LatticeModel.h"
using
namespace
walberla
;
typedef
lbm
::
FluctuatingMRT_LatticeModel
LatticeModel_T
;
typedef
LatticeModel_T
::
Stencil
Stencil_T
;
typedef
LatticeModel_T
::
CommunicationStencil
CommunicationStencil_T
;
typedef
lbm
::
PdfField
<
LatticeModel_T
>
PdfField_T
;
typedef
GhostLayerField
<
real_t
,
LatticeModel_T
::
Stencil
::
D
>
VectorField_T
;
typedef
GhostLayerField
<
real_t
,
1
>
ScalarField_T
;
typedef
walberla
::
uint8_t
flag_t
;
typedef
FlagField
<
flag_t
>
FlagField_T
;
int
main
(
int
argc
,
char
**
argv
)
{
walberla
::
Environment
walberlaEnv
(
argc
,
argv
);
auto
blocks
=
blockforest
::
createUniformBlockGridFromConfig
(
walberlaEnv
.
config
()
);
// read parameters
auto
parameters
=
walberlaEnv
.
config
()
->
getOneBlock
(
"Parameters"
);
const
real_t
omega
=
parameters
.
getParameter
<
real_t
>
(
"omega"
,
real_c
(
1.4
)
);
const
Vector3
<
real_t
>
initialVelocity
=
parameters
.
getParameter
<
Vector3
<
real_t
>
>
(
"initialVelocity"
,
Vector3
<
real_t
>
()
);
const
uint_t
timesteps
=
parameters
.
getParameter
<
uint_t
>
(
"timesteps"
,
uint_c
(
10
)
);
const
real_t
temperature
=
real_t
(
0.01
);
const
uint_t
seed
=
uint_t
(
0
);
const
double
remainingTimeLoggerFrequency
=
parameters
.
getParameter
<
double
>
(
"remainingTimeLoggerFrequency"
,
3.0
);
// in seconds
// create fields
BlockDataID
forceFieldId
=
field
::
addToStorage
<
VectorField_T
>
(
blocks
,
"Force"
,
real_t
(
0.0
));
LatticeModel_T
latticeModel
=
LatticeModel_T
(
forceFieldId
,
omega
,
omega
,
seed
,
temperature
,
uint_t
(
0
)
);
BlockDataID
pdfFieldId
=
lbm
::
addPdfFieldToStorage
(
blocks
,
"pdf field"
,
latticeModel
,
initialVelocity
,
real_t
(
1
)
);
BlockDataID
flagFieldId
=
field
::
addFlagFieldToStorage
<
FlagField_T
>
(
blocks
,
"flag field"
);
// create and initialize boundary handling
const
FlagUID
fluidFlagUID
(
"Fluid"
);
// create time loop
SweepTimeloop
timeloop
(
blocks
->
getBlockStorage
(),
timesteps
);
// create communication for PdfField
blockforest
::
communication
::
UniformBufferedScheme
<
CommunicationStencil_T
>
communication
(
blocks
);
communication
.
addPackInfo
(
make_shared
<
lbm
::
PdfFieldPackInfo
<
LatticeModel_T
>
>
(
pdfFieldId
)
);
// add LBM sweep and communication to time loop
timeloop
.
add
()
<<
BeforeFunction
(
[
&
](){
latticeModel
.
time_step
=
uint32_c
(
timeloop
.
getCurrentTimeStep
());
},
"set RNG counter"
)
<<
BeforeFunction
(
communication
,
"communication"
)
<<
Sweep
(
LatticeModel_T
::
Sweep
(
pdfFieldId
),
"LB stream & collide"
);
// LBM stability check
timeloop
.
addFuncAfterTimeStep
(
makeSharedFunctor
(
field
::
makeStabilityChecker
<
PdfField_T
,
FlagField_T
>
(
walberlaEnv
.
config
(),
blocks
,
pdfFieldId
,
flagFieldId
,
fluidFlagUID
)
),
"LBM stability check"
);
// log remaining time
timeloop
.
addFuncAfterTimeStep
(
timing
::
RemainingTimeLogger
(
timeloop
.
getNrOfTimeSteps
(),
remainingTimeLoggerFrequency
),
"remaining time logger"
);
// add VTK output to time loop
lbm
::
VTKOutput
<
LatticeModel_T
,
FlagField_T
>::
addToTimeloop
(
timeloop
,
blocks
,
walberlaEnv
.
config
(),
pdfFieldId
,
flagFieldId
,
fluidFlagUID
);
// create adaptors, so that the GUI also displays density and velocity
// adaptors are like fields with the difference that they do not store values
// but calculate the values based on other fields ( here the PdfField )
field
::
addFieldAdaptor
<
lbm
::
Adaptor
<
LatticeModel_T
>::
Density
>
(
blocks
,
pdfFieldId
,
"DensityAdaptor"
);
field
::
addFieldAdaptor
<
lbm
::
Adaptor
<
LatticeModel_T
>::
VelocityVector
>
(
blocks
,
pdfFieldId
,
"VelocityAdaptor"
);
if
(
parameters
.
getParameter
<
bool
>
(
"useGui"
,
false
)
)
{
GUI
gui
(
timeloop
,
blocks
,
argc
,
argv
);
lbm
::
connectToGui
<
LatticeModel_T
>
(
gui
);
gui
.
run
();
}
else
timeloop
.
run
();
return
EXIT_SUCCESS
;
}
tests/lbm/codegen/FluctuatingMRT.py
0 → 100644
View file @
eb1fd933
import
sympy
as
sp
import
pystencils
as
ps
from
lbmpy.creationfunctions
import
create_lb_collision_rule
from
pystencils_walberla
import
CodeGeneration
from
lbmpy_walberla
import
generate_lattice_model
with
CodeGeneration
()
as
ctx
:
omega_shear
,
omega_bulk
=
sp
.
symbols
(
"omega_shear, omega_bulk"
)
temperature
=
sp
.
symbols
(
"temperature"
)
force_field
,
vel_field
=
ps
.
fields
(
"force(3), velocity(3): [3D]"
,
layout
=
'fzyx'
)
collision_rule
=
create_lb_collision_rule
(
stencil
=
'D3Q19'
,
compressible
=
True
,
fluctuating
=
{
'temperature'
:
temperature
,
'block_offsets'
:
'walberla'
,
},
method
=
'mrt3'
,
relaxation_rates
=
[
omega_shear
,
omega_bulk
],
force_model
=
'guo'
,
force
=
force_field
.
center_vector
,
optimization
=
{
'cse_global'
:
True
}
)
generate_lattice_model
(
ctx
,
'FluctuatingMRT_LatticeModel'
,
collision_rule
)
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment