OpenLB 1.7
Loading...
Searching...
No Matches
setFreeEnergyInletBoundary2D.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 Free Energy Inlet Boundary
25//This is a new version of the Boundary, which only contains free floating functions
26#ifndef SET_FREE_ENERGY_INLET_BOUNDARY_2D_HH
27#define SET_FREE_ENERGY_INLET_BOUNDARY_2D_HH
28
30
31namespace olb {
32
34template<typename T, typename DESCRIPTOR, typename MixinDynamics>
35void setFreeEnergyInletBoundary(SuperLattice<T, DESCRIPTOR>& sLattice, T omega, SuperGeometry<T,2>& 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<SuperIndicatorF2D<T>>&& indicator, std::string type, int latticeNumber)
43{
44 OstreamManager clout(std::cout, "setFreeEnergyInletBoundary");
45 bool includeOuterCells = false;
46 int _overlap = 1;
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<T,DESCRIPTOR>(sLattice, std::forward<decltype(indicator)>(indicator), _overlap);
57}
58
60template<typename T, typename DESCRIPTOR, typename MixinDynamics>
61void setFreeEnergyInletBoundary(BlockLattice<T,DESCRIPTOR>& block,T omega, BlockIndicatorF2D<T>& indicator, std::string type,
62 int latticeNumber, bool includeOuterCells)
63{
64 bool _output = false;
65 OstreamManager clout(std::cout, "setFreeEnergyInletBoundary");
66 auto& blockGeometryStructure = indicator.getBlockGeometry();
67 const int margin = includeOuterCells ? 0 : 1;
68 std::vector<int> discreteNormal(3, 0);
69 blockGeometryStructure.forSpatialLocations([&](auto iX, auto iY) {
70 if (blockGeometryStructure.getNeighborhoodRadius({iX, iY}) >= margin
71 && indicator(iX, iY)) {
72 discreteNormal = blockGeometryStructure.getStatistics().getType(iX,iY);
73 if (discreteNormal[0] == 0) {
74 Dynamics<T, DESCRIPTOR>* dynamics = nullptr;
75 if (discreteNormal[1] == -1) {
76
77 if (latticeNumber == 1) {
78 //set momenta and dynamics for a pressure boundary on indicated cells
79 if (type == "density") {
80 dynamics = block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
81 MixinDynamics, momenta::RegularizedPressureBoundaryTuple<0,-1>>>();
82 }
83 else {
84 //set momenta and dynamics for a velocity boundary on indicated cells
85 dynamics = block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
86 MixinDynamics, momenta::RegularizedVelocityBoundaryTuple<0,-1>>>();
87 }
88 block.defineDynamics({iX, iY}, dynamics);
89 }
90 else {
91 dynamics = block.template getDynamics<FreeEnergyInletOutletDynamics<T,DESCRIPTOR,0,-1>>();
92 block.defineDynamics({iX, iY}, dynamics);
93 }
94 }
95
96 else if (discreteNormal[1] == 1) {
97 if (latticeNumber == 1) {
98 if (type == "density") {
99 //set momenta and dynamics for a pressure boundary on indicated cells
100 dynamics = block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
102 }
103 else {
104 //set momenta and dynamics for a velocity boundary on indicated cells
105 dynamics = block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
107 }
108 block.defineDynamics({iX, iY}, dynamics);
109 }
110 else {
111 dynamics = block.template getDynamics<FreeEnergyInletOutletDynamics<T,DESCRIPTOR,0,1>>();
112 block.defineDynamics({iX, iY}, dynamics);
113 }
114 }
115
116 else if (discreteNormal[2] == -1) {
117 if (latticeNumber == 1) {
118 if (type == "density") {
119 //set momenta and dynamics for a pressure boundary on indicated cells
120 dynamics = block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
121 MixinDynamics, momenta::RegularizedPressureBoundaryTuple< 1,-1>>>();
122 }
123 else {
124 //set momenta and dynamics for a velocity boundary on indicated cells
125 dynamics = block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR,
126 MixinDynamics, momenta::RegularizedVelocityBoundaryTuple< 1,-1>>>();
127 }
128 block.defineDynamics({iX, iY}, dynamics);
129 }
130 else {
131 dynamics = block.template getDynamics<FreeEnergyInletOutletDynamics<T,DESCRIPTOR,1,-1>>();
132 block.defineDynamics({iX, iY}, dynamics);
133 }
134 }
135
136 else if (discreteNormal[2] == 1) {
137 if (latticeNumber == 1) {
138 if (type == "density") {
139 //set momenta and dynamics for a pressure boundary on indicated cells
140 dynamics = block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR, MixinDynamics,
142 }
143 else {
144 //set momenta and dynamics for a velocity boundary on indicated cells
145 dynamics = block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR, MixinDynamics,
147 }
148 block.defineDynamics({iX, iY}, dynamics);
149 }
150 else {
151 dynamics = block.template getDynamics<FreeEnergyInletOutletDynamics<T,DESCRIPTOR,1,1>>();
152 block.defineDynamics({iX, iY}, dynamics);
153 }
154 }
155
156 if (latticeNumber != 1) {
157 block.defineDynamics({iX, iY}, dynamics);
158 }
159
160 if ( latticeNumber == 1){
161 block.addPostProcessor(
162 typeid(stage::PostStream), {iX, iY},
163 olb::boundaryhelper::promisePostProcessorForNormal<T,DESCRIPTOR, FreeEnergyChemPotBoundaryProcessor2DA>(
164 Vector<int,2>(discreteNormal.data() + 1)));
165 } else {
166 block.addPostProcessor(
167 typeid(stage::PostStream), {iX, iY},
168 olb::boundaryhelper::promisePostProcessorForNormal<T,DESCRIPTOR, FreeEnergyChemPotBoundaryProcessor2DB>(
169 Vector<int,2>(discreteNormal.data() + 1)));
170 }
171
172 dynamics->getParameters(block).template set<descriptors::OMEGA>(omega);
173 //sets the boundary on any indicated cell
174 //located in setLocalVelocityBoundary2D.h/hh
175 setBoundary(block, iX,iY, dynamics);
176
177 if (_output) {
178 clout << "setFreeEnergyInletBoundary<" << "," << ">(" << iX << ", "<< iY << ")" << std::endl;
179 }
180
181 }
182 }
183 });
184}
185
186}//namespace olb
187#endif
Base block indicator functor (discrete)
BlockGeometry< T, 2 > & 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.
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.
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 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...
void setBoundary(BlockLattice< T, DESCRIPTOR > &block, int iX, int iY, Dynamics< T, DESCRIPTOR > *dynamics, PostProcessorGenerator2D< T, DESCRIPTOR > *postProcessor)
Interface for per-cell dynamics.
Definition interface.h:54
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