OpenLB 1.7
Loading...
Searching...
No Matches
setFreeEnergyInletBoundary3D.hh
Go to the documentation of this file.
1/* This file is part of the OpenLB library
2 *
3 * Copyright (C) 2020 Alexander Schulz
4 * E-mail contact: info@openlb.net
5 * The most recent release of OpenLB can be downloaded at
6 * <http://www.openlb.net/>
7 *
8 * This program is free software; you can redistribute it and/or
9 * modify it under the terms of the GNU General Public License
10 * as published by the Free Software Foundation; either version 2
11 * of the License, or (at your option) any later version.
12 *
13 * This program is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public
19 * License along with this program; if not, write to the Free
20 * Software Foundation, Inc., 51 Franklin Street, Fifth Floor,
21 * Boston, MA 02110-1301, USA.
22*/
23
24//This file contains the FreeEnergyInlet Boundary
25//This is a new version of the Boundary, which only contains free floating functions
26#ifndef SET_FREE_ENERGY_INLET_BOUNDARY_3D_HH
27#define SET_FREE_ENERGY_INLET_BOUNDARY_3D_HH
28
30
31namespace olb {
32
34template<typename T, typename DESCRIPTOR, typename MixinDynamics>
35void setFreeEnergyInletBoundary(SuperLattice<T, DESCRIPTOR>& sLattice,T omega, SuperGeometry<T,3>& superGeometry, int material, std::string type, int latticeNumber)
36{
37 setFreeEnergyInletBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice, omega, superGeometry.getMaterialIndicator(material), type, latticeNumber);
38}
39
41template<typename T, typename DESCRIPTOR, typename MixinDynamics>
42void setFreeEnergyInletBoundary(SuperLattice<T, DESCRIPTOR>& sLattice,T omega,FunctorPtr<SuperIndicatorF3D<T>>&& indicator,std::string type, int latticeNumber)
43{
44 bool includeOuterCells = false;
45 int _overlap = 1;
46 OstreamManager clout(std::cout, "setFreeEnergyInletBoundary");
47 if (indicator->getSuperGeometry().getOverlap() == 1) {
48 includeOuterCells = true;
49 clout << "WARNING: overlap == 1, boundary conditions set on overlap despite unknown neighbor materials" << std::endl;
50 }
51 for (int iCloc = 0; iCloc < sLattice.getLoadBalancer().size(); ++iCloc) {
52 setFreeEnergyInletBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice.getBlock(iCloc),omega,
53 indicator->getBlockIndicatorF(iCloc), type, latticeNumber, includeOuterCells);
54 }
56 addPoints2CommBC(sLattice, std::forward<decltype(indicator)>(indicator), _overlap);
57}
58
59
61template<typename T, typename DESCRIPTOR, typename MixinDynamics>
62void setFreeEnergyInletBoundary(BlockLattice<T,DESCRIPTOR>& _block, T omega, BlockIndicatorF3D<T>& indicator, std::string type,
63 int latticeNumber, bool includeOuterCells)
64{
65 auto& blockGeometryStructure = indicator.getBlockGeometry();
66 const int margin = includeOuterCells ? 0 : 1;
67 OstreamManager clout(std::cout, "setFreeEnergyInletBoundary");
68 bool _output = false;
69 std::vector<int> discreteNormal(4, 0);
70 blockGeometryStructure.forSpatialLocations([&](auto iX, auto iY, auto iZ) {
71 if (blockGeometryStructure.getNeighborhoodRadius({iX, iY, iZ}) >= margin
72 && indicator(iX, iY, iZ)) {
73 discreteNormal = blockGeometryStructure.getStatistics().getType(iX, iY, iZ, true);
74
75 if (discreteNormal[0] == 0) {
76 Dynamics<T, DESCRIPTOR>* dynamics = nullptr;
77 if (discreteNormal[1] != 0 && discreteNormal[1] == -1) {
78 if (latticeNumber == 1) {//set PressureBoundary momenta and dynamics on indicated cells
79 if (type == "density") {
80 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
81 MixinDynamics,momenta::RegularizedPressureBoundaryTuple<0,-1>>>();
82 }
83 else { //set VelocityBoundary momenta and dynamics on indicated cells
84 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
86 >>();
87 }
88 }
89 else {
90 dynamics = _block.template getDynamics<FreeEnergyInletOutletDynamics<T,DESCRIPTOR,0,-1>>();
91 }
92 }
93
94 else if (discreteNormal[1] != 0 && discreteNormal[1] == 1) {
95 if (latticeNumber == 1) {
96 if (type == "density") {//set PressureBoundary momenta and dynamics on indicated cells
97 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
99 }
100 else { //set VelocityBoundary momenta and dynamics on indicated cells
101 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
103 >>();
104 }
105 }
106 else {
107 dynamics = _block.template getDynamics<FreeEnergyInletOutletDynamics<T,DESCRIPTOR,0,1>>();
108 }
109 }
110
111 else if (discreteNormal[2] != 0 && discreteNormal[2] == -1) {
112 if (latticeNumber == 1) {
113 if (type == "density") {//set PressureBoundary momenta and dynamics on indicated cells
114 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
115 MixinDynamics,momenta::RegularizedPressureBoundaryTuple<1,-1>>>();
116 }
117 else { //set VelocityBoundary momenta and dynamics on indicated cells
118 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
120 >>();
121 }
122 }
123 else {
124 dynamics = _block.template getDynamics<FreeEnergyInletOutletDynamics<T,DESCRIPTOR,1,-1>>();
125 }
126 }
127
128 else if (discreteNormal[2] != 0 && discreteNormal[2] == 1) {
129 if (latticeNumber == 1) {
130 if (type == "density") {//set PressureBoundary momenta and dynamics on indicated cells
131 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
133 }
134 else { //set VelocityBoundary momenta and dynamics on indicated cells
135 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
137 >>();
138 }
139 }
140 else {
141 dynamics = _block.template getDynamics<FreeEnergyInletOutletDynamics<T,DESCRIPTOR,1,1>>();
142 }
143 }
144
145 else if (discreteNormal[3] != 0 && discreteNormal[3] == -1) {
146 if (latticeNumber == 1) {
147 if (type == "density") {//set PressureBoundary momenta and dynamics on indicated cells
148 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
149 MixinDynamics,momenta::RegularizedPressureBoundaryTuple<2,-1>>>();
150 }
151 else { //set VelocityBoundary momenta and dynamics on indicated cells
152 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
154 >>();
155 }
156 }
157 else {
158 dynamics = _block.template getDynamics<FreeEnergyInletOutletDynamics<T,DESCRIPTOR,2,-1>>();
159 }
160 }
161
162 else if (discreteNormal[3] != 0 && discreteNormal[3] == 1) {
163 if (latticeNumber == 1) {
164 if (type == "density") {//set PressureBoundary momenta and dynamics on indicated cells
165 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
167 }
168 else { //set VelocityBoundary momenta and dynamics on indicated cells
169 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
171 >>();
172 }
173 }
174 else {
175 dynamics = _block.template getDynamics<FreeEnergyInletOutletDynamics<T,DESCRIPTOR,2,1>>();
176 }
177 }
178 dynamics->getParameters(_block).template set<descriptors::OMEGA>(omega);
179 if (latticeNumber != 1) {
180 _block.defineDynamics({iX, iY, iZ}, dynamics);
181 auto cell = _block.get(iX,iY,iZ);
182 dynamics->initialize(cell);
183 }
184 }
185
186 if (latticeNumber == 1){
187 _block.addPostProcessor(
188 typeid(stage::PostStream), {iX, iY, iZ},
189 olb::boundaryhelper::promisePostProcessorForNormal<T,DESCRIPTOR, FreeEnergyChemPotBoundaryProcessor3DA>(
190 Vector<int,3>(discreteNormal.data() + 1)));
191 } else {
192 _block.addPostProcessor(
193 typeid(stage::PostStream), {iX, iY, iZ},
194 olb::boundaryhelper::promisePostProcessorForNormal<T,DESCRIPTOR, FreeEnergyChemPotBoundaryProcessor3DB>(
195 Vector<int,3>(discreteNormal.data() + 1)));
196 }
197
198 if (_output) {
199 clout << "setFreeEnergyInletBoundary<" << "," << ">(" << iX << ", "<< iY << ", " << iZ << ")" << std::endl;
200 }
201 }
202 });
203}
204
205
206
207}//namespace olb
208#endif
Base block indicator functor.
BlockGeometry< T, 3 > & getBlockGeometry()
Get underlying block geometry structure.
Platform-abstracted block lattice for external access and inter-block interaction.
virtual void addPostProcessor(std::type_index stage, LatticeR< DESCRIPTOR::d > latticeR, PostProcessorPromise< T, DESCRIPTOR > &&promise)=0
Schedule post processor for application to latticeR in stage.
void defineDynamics(LatticeR< DESCRIPTOR::d > latticeR, DynamicsPromise< T, DESCRIPTOR > &&promise)
Assign promised DYNAMICS to latticeR.
Cell< T, DESCRIPTOR > get(CellID iCell)
Get Cell interface for index iCell.
Regularized BGK collision followed by any other Dynamics.
Definition dynamics.h:193
Smart pointer for managing the various ways of passing functors around.
Definition functorPtr.h:60
class for marking output with some text
Representation of a statistic for a parallel 2D geometry.
std::unique_ptr< SuperIndicatorF< T, D > > getMaterialIndicator(std::vector< int > &&materials)
Returns a material indicator using the given vector of materials.
Base indicator functor (discrete)
Super class maintaining block lattices for a cuboid decomposition.
BlockLattice< T, DESCRIPTOR > & getBlock(int locC)
Return BlockLattice with local index locC.
LoadBalancer< T > & getLoadBalancer()
Read and write access to the load balancer.
Plain old scalar vector.
Definition vector.h:47
Top level namespace for all of OpenLB.
void addPoints2CommBC(SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF2D< T > > &&indicator, int _overlap)
Adds needed Cells to the Communicator _commBC in SuperLattice.
void setFreeEnergyInletBoundary(SuperLattice< T, DESCRIPTOR > &sLattice, T omega, SuperGeometry< T, 2 > &superGeometry, int material, std::string type, int latticeNumber)
Implementation of a inlet boundary condition for the partner lattices of the binary or ternary free e...
Interface for per-cell dynamics.
Definition interface.h:54
virtual void initialize(Cell< T, DESCRIPTOR > &cell)
Initialize dynamics-specific data for cell.
Definition interface.h:68
virtual AbstractParameters< T, DESCRIPTOR > & getParameters(BlockLattice< T, DESCRIPTOR > &block)=0
Parameters access for legacy post processors.
Dynamics constructed as a tuple of momenta, equilibrium and collision.
Definition interface.h:182
Communication after propagation.
Definition stages.h:36