OpenLB 1.7
Loading...
Searching...
No Matches
setLocalVelocityBoundary3D.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 Local Velocity Boundary
25//This is a new version of the Boundary, which only contains free floating functions
26#ifndef SET_LOCAL_VELOCITY_BOUNDARY_HH
27#define SET_LOCAL_VELOCITY_BOUNDARY_HH
28
30
31namespace olb {
32
33
35template<typename T,typename DESCRIPTOR, class MixinDynamics>
36void setLocalVelocityBoundary(SuperLattice<T, DESCRIPTOR>& sLattice, T omega, SuperGeometry<T,3>& superGeometry, int material)
37{
38 setLocalVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice, omega, superGeometry.getMaterialIndicator(material));
39}
40
42template<typename T, typename DESCRIPTOR, class MixinDynamics>
44{
45 OstreamManager clout(std::cout, "setLocalVelocityBoundary");
46 int _overlap = 1;
47 bool includeOuterCells = false;
48 if (indicator->getSuperGeometry().getOverlap() == 1) {
49 includeOuterCells = true;
50 clout << "WARNING: overlap == 1, boundary conditions set on overlap despite unknown neighbor materials" << std::endl;
51 }
52
53 for (int iC = 0; iC < sLattice.getLoadBalancer().size(); ++iC) {
54 setLocalVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice.getBlock(iC), omega, indicator->getBlockIndicatorF(iC),includeOuterCells);
55 }
57 addPoints2CommBC(sLattice, std::forward<decltype(indicator)>(indicator), _overlap);
58}
59
60
61
63template<typename T, typename DESCRIPTOR, class MixinDynamics>
64void setLocalVelocityBoundary(BlockLattice<T,DESCRIPTOR>& _block,T omega, BlockIndicatorF3D<T>& indicator, bool includeOuterCells)
65{
66 using namespace boundaryhelper;
67 const auto& blockGeometryStructure = indicator.getBlockGeometry();
68 const int margin = includeOuterCells ? 0 : 1;
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 Dynamics<T,DESCRIPTOR>* dynamics = nullptr;
74 discreteNormal = indicator.getBlockGeometry().getStatistics().getType(iX, iY, iZ);
75 if (discreteNormal[0] == 0) {
76 //setVelocityBoundary
77 if (discreteNormal[1] != 0 && discreteNormal[1] == -1) {
78 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR, MixinDynamics,
80 >>();
81 }
82 else if (discreteNormal[1] != 0 && discreteNormal[1] == 1) {
83 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR, MixinDynamics,
85 >>();
86 }
87 else if (discreteNormal[2] != 0 && discreteNormal[2] == -1) {
88 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR, MixinDynamics,
90 >>();
91 }
92 else if (discreteNormal[2] != 0 && discreteNormal[2] == 1) {
93 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR, MixinDynamics,
95 >>();
96 }
97 else if (discreteNormal[3] != 0 && discreteNormal[3] == -1) {
98 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR, MixinDynamics,
100 >>();
101 }
102 else if (discreteNormal[3] != 0 && discreteNormal[3] == 1) {
103 dynamics = _block.template getDynamics<CombinedRLBdynamics<T,DESCRIPTOR, MixinDynamics,
105 >>();
106 }
107 }
108 //ExternalVelocityCorner
109 else if (discreteNormal[0] == 1) {
110 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
112 >>();
113 _block.addPostProcessor(
114 typeid(stage::PostStream), {iX, iY, iZ},
115 promisePostProcessorForNormal<T,DESCRIPTOR,OuterVelocityCornerProcessor3D>(
116 Vector<int,3>(discreteNormal.data() + 1)));
117 }
118 //InternalVelocityCorner
119 else if (discreteNormal[0] == 2) {
120 dynamics = _block.getDynamics(PlainMixinDynamicsForNormalMomenta<T,DESCRIPTOR,
122 >::construct(Vector<int,3>(discreteNormal.data() + 1)));
123 }
124 //ExternalVelocityEdge
125 else if (discreteNormal[0] == 3) {
126 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
128 >>();
129 _block.addPostProcessor(
130 typeid(stage::PostStream), {iX, iY, iZ},
131 promisePostProcessorForNormalSpecial<T,DESCRIPTOR,OuterVelocityEdgeProcessor3D>(
132 Vector<int,3>(discreteNormal.data() + 1)));
133 }
134 //InternalVelocityEdge
135 else if (discreteNormal[0] == 4) {
136 dynamics = _block.getDynamics(PlainMixinDynamicsForNormalSpecialMomenta<T,DESCRIPTOR,
138 >::construct(Vector<int,3>(discreteNormal.data() + 1)));
139 }
140 if (dynamics) {
141 dynamics->getParameters(_block).template set<descriptors::OMEGA>(omega);
142 }
143 setBoundary(_block, iX,iY,iZ, dynamics);
144 }
145 });
146}
147
148}//namespace olb
149
150#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.
virtual Dynamics< T, DESCRIPTOR > * getDynamics(DynamicsPromise< T, DESCRIPTOR > &&)=0
Return pointer to dynamics yielded by promise.
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 setLocalVelocityBoundary(SuperLattice< T, DESCRIPTOR > &sLattice, T omega, SuperGeometry< T, 2 > &superGeometry, int material)
Initialising the setLocalVelocityBoundary function on the superLattice domain.
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.
Communication after propagation.
Definition stages.h:36