OpenLB 1.7
Loading...
Searching...
No Matches
setInterpolatedVelocityBoundary3D.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 Interpolated Velocity Boundary
25//This is a new version of the Boundary, which only contains free floating functions
26#ifndef SET_INTERPOLATED_VELOCITY_BOUNDARY_HH
27#define SET_INTERPOLATED_VELOCITY_BOUNDARY_HH
28
30
31namespace olb {
32
34template<typename T,typename DESCRIPTOR, class MixinDynamics>
35void setInterpolatedVelocityBoundary(SuperLattice<T, DESCRIPTOR>& sLattice, T omega, SuperGeometry<T,3>& superGeometry, int material)
36{
37 setInterpolatedVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice, omega, superGeometry.getMaterialIndicator(material));
38}
39
41template<typename T, typename DESCRIPTOR, class MixinDynamics>
43{
44 OstreamManager clout(std::cout, "setInterpolatedVelocityBoundary");
45 int _overlap = 1;
46 bool includeOuterCells = false;
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 iC = 0; iC < sLattice.getLoadBalancer().size(); ++iC) {
52 setInterpolatedVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice.getBlock(iC),omega, indicator->getBlockIndicatorF(iC),includeOuterCells);
53 }
55 //the addPoints2CommBC function is initialised inside setLocalVelocityBoundary3D.h/hh
56 addPoints2CommBC(sLattice,std::forward<decltype(indicator)>(indicator), _overlap);
57}
58
59
60//set InterpolatedVelocityBoundary on indicated cells inside the block domain
61template<typename T, typename DESCRIPTOR, class MixinDynamics>
62void setInterpolatedVelocityBoundary(BlockLattice<T,DESCRIPTOR>& _block, T omega,BlockIndicatorF3D<T>& indicator, bool includeOuterCells)
63{
64 using namespace boundaryhelper;
65 auto& blockGeometryStructure = indicator.getBlockGeometry();
66 const int margin = includeOuterCells ? 0 : 1;
67 std::vector<int> discreteNormal(4,0);
68 blockGeometryStructure.forSpatialLocations([&](auto iX, auto iY, auto iZ) {
69 if (blockGeometryStructure.getNeighborhoodRadius({iX, iY, iZ}) >= margin
70 && indicator(iX, iY, iZ)) {
71 Dynamics<T,DESCRIPTOR>* dynamics = nullptr;
72 discreteNormal = indicator.getBlockGeometry().getStatistics().getType(iX, iY, iZ);
73 if (discreteNormal[0] == 0) {
74 if (discreteNormal[1] != 0 && discreteNormal[1] == -1) {//set momenta, dynamics and postProcessors on indicated velocityBoundaryCells
75 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
77 >>();
78 _block.addPostProcessor(
79 typeid(stage::PostStream), {iX, iY, iZ},
80 meta::id<PlaneFdBoundaryProcessor3D<T,DESCRIPTOR,0,-1>>());
81 }
82 else if (discreteNormal[1] != 0 && discreteNormal[1] == 1) {
83 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
85 >>();
86 _block.addPostProcessor(
87 typeid(stage::PostStream), {iX, iY, iZ},
89 }
90 else if (discreteNormal[2] != 0 && discreteNormal[2] == -1) {
91 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
93 >>();
94 _block.addPostProcessor(
95 typeid(stage::PostStream), {iX, iY, iZ},
96 meta::id<PlaneFdBoundaryProcessor3D<T,DESCRIPTOR,1,-1>>());
97 }
98 else if (discreteNormal[2] != 0 && discreteNormal[2] == 1) {
99 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
101 >>();
102 _block.addPostProcessor(
103 typeid(stage::PostStream), {iX, iY, iZ},
105 }
106 else if (discreteNormal[3] != 0 && discreteNormal[3] == -1) {
107 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
109 >>();
110 _block.addPostProcessor(
111 typeid(stage::PostStream), {iX, iY, iZ},
112 meta::id<PlaneFdBoundaryProcessor3D<T,DESCRIPTOR,2,-1>>());
113 }
114 else if (discreteNormal[3] != 0 && discreteNormal[3] == 1) {
115 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
117 >>();
118 _block.addPostProcessor(
119 typeid(stage::PostStream), {iX, iY, iZ},
121 }
122 }
123
124 else if (discreteNormal[0] == 1) {//set momenta,dynamics and postProcessors on indicated velocityBoundary External Corner cells
125 //ExternalVelocityCorner
126 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
128 >>();
129 _block.addPostProcessor(
130 typeid(stage::PostStream), {iX, iY, iZ},
131 promisePostProcessorForNormal<T,DESCRIPTOR,OuterVelocityCornerProcessor3D>(
132 Vector<int,3>(discreteNormal.data() + 1)));
133 }
134
135 else if (discreteNormal[0] == 2) {//
136 //Internalvelocitycorner
137 dynamics = _block.getDynamics(PlainMixinDynamicsForNormalMomenta<T,DESCRIPTOR,
139 >::construct(Vector<int,3>(discreteNormal.data() + 1)));
140 }
141 //ExternalVelocityEdge
142 else if (discreteNormal[0] == 3) {//set momenta,dynamics and postProcessors on indicated velocityBoundary External Edge cells
143 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
145 >>();
146 _block.addPostProcessor(
147 typeid(stage::PostStream), {iX, iY, iZ},
148 promisePostProcessorForNormalSpecial<T,DESCRIPTOR,OuterVelocityEdgeProcessor3D>(
149 Vector<int,3>(discreteNormal.data() + 1)));
150 }
151 //InternalVelocityEdge
152 else if (discreteNormal[0] == 4) {//set momenta,dynamics and postProcessors on indicated velocityBoundary Inner Edge cells
153 dynamics = _block.getDynamics(PlainMixinDynamicsForNormalSpecialMomenta<T,DESCRIPTOR,
155 >::construct(Vector<int,3>(discreteNormal.data() + 1)));
156 }
157
158 if (dynamics) {
159 dynamics->getParameters(_block).template set<descriptors::OMEGA>(omega);
160 }
161 setBoundary(_block, iX,iY,iZ, dynamics);
162 }
163 });
164}
165
166}
167#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
This class computes the skordos BC on a plane wall in 3D but with a limited number of terms added to ...
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 setInterpolatedVelocityBoundary(SuperLattice< T, DESCRIPTOR > &sLattice, T omega, SuperGeometry< T, 2 > &superGeometry, int material)
Initialising the setInterpolatedVelocityBoundary function on the superLattice domain Interpolated Bou...
void addPoints2CommBC(SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF2D< T > > &&indicator, int _overlap)
Adds needed Cells to the Communicator _commBC in SuperLattice.
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.
Identity type to pass non-constructible types as value.
Definition meta.h:79
Communication after propagation.
Definition stages.h:36