OpenLB 1.7
Loading...
Searching...
No Matches
setInterpolatedPressureBoundary3D.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 Pressure Boundary
25//This is an onLattice Boundary
26//This is a new version of the Boundary, which only contains free floating functions
27#ifndef SET_INTERPOLATED_PRESSURE_BOUNDARY_HH
28#define SET_INTERPOLATED_PRESSURE_BOUNDARY_HH
29
31
32namespace olb {
33
35template<typename T,typename DESCRIPTOR, class MixinDynamics>
37 int material)
38{
39 setInterpolatedPressureBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice, omega, superGeometry.getMaterialIndicator(material));
40}
41
43template<typename T, typename DESCRIPTOR, typename MixinDynamics>
45{
46 OstreamManager clout(std::cout, "setInterpolatedPressureBoundary");
47 int _overlap = 1;
48 bool includeOuterCells = false;
49 if (indicator->getSuperGeometry().getOverlap() == 1) {
50 includeOuterCells = true;
51 clout << "WARNING: overlap == 1, boundary conditions set on overlap despite unknown neighbor materials" << std::endl;
52 }
53 //clout << sLattice.getLoadBalancer().size() <<"sLattice.getLoadBalancer.size()" << std::endl;
54 for (int iC = 0; iC < sLattice.getLoadBalancer().size(); ++iC) {
55 setInterpolatedPressureBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice.getBlock(iC), omega, indicator->getBlockIndicatorF(iC),
56 includeOuterCells);
57 }
58 //defined in setLocalVelocityBoundary
59 addPoints2CommBC(sLattice, std::forward<decltype(indicator)>(indicator), _overlap);
60
61}
62
63
65template<typename T, typename DESCRIPTOR, typename MixinDynamics>
66void setInterpolatedPressureBoundary(BlockLattice<T,DESCRIPTOR>& _block, T omega, BlockIndicatorF3D<T>& indicator, bool includeOuterCells)
67{
68 auto& blockGeometryStructure = indicator.getBlockGeometry();
69 const int margin = includeOuterCells ? 0 : 1;
70 std::vector<int> discreteNormal(4,0);
71 blockGeometryStructure.forSpatialLocations([&](auto iX, auto iY, auto iZ) {
72 if (blockGeometryStructure.getNeighborhoodRadius({iX, iY, iZ}) >= margin
73 && indicator(iX, iY, iZ)) {
74 Dynamics<T,DESCRIPTOR>* dynamics = nullptr;
75 discreteNormal = blockGeometryStructure.getStatistics().getType(iX, iY, iZ);
76
77 if (discreteNormal[0] == 0) {
78 if (discreteNormal[1] != 0 && discreteNormal[1] == -1) {
79 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
81 >>();
82 _block.addPostProcessor(
83 typeid(stage::PostStream), {iX, iY, iZ},
84 meta::id<PlaneFdBoundaryProcessor3D<T,DESCRIPTOR,0,-1>>());
85 }
86 else if (discreteNormal[1] != 0 && discreteNormal[1] == 1) {
87 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
89 >>();
90 _block.addPostProcessor(
91 typeid(stage::PostStream), {iX, iY, iZ},
93 }
94 else if (discreteNormal[2] != 0 && discreteNormal[2] == -1) {
95 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
97 >>();
98 _block.addPostProcessor(
99 typeid(stage::PostStream), {iX, iY, iZ},
100 meta::id<PlaneFdBoundaryProcessor3D<T,DESCRIPTOR,1,-1>>());
101 }
102 else if (discreteNormal[2] != 0 && discreteNormal[2] == 1) {
103 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
105 >>();
106 _block.addPostProcessor(
107 typeid(stage::PostStream), {iX, iY, iZ},
109 }
110 else if (discreteNormal[3] != 0 && discreteNormal[3] == -1) {
111 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
113 >>();
114 _block.addPostProcessor(
115 typeid(stage::PostStream), {iX, iY, iZ},
116 meta::id<PlaneFdBoundaryProcessor3D<T,DESCRIPTOR,2,-1>>());
117 }
118 else if (discreteNormal[3] != 0 && discreteNormal[3] == 1) {
119 dynamics = _block.template getDynamics<typename MixinDynamics::template exchange_momenta<
121 >>();
122 _block.addPostProcessor(
123 typeid(stage::PostStream), {iX, iY, iZ},
125 }
126 else{
127 std::cout << "Discrete Normal is not valid for interpolated pressure boundary!" << std::endl;
129 }
130 }
131 else {
132 // When this error occurs, check your geometry for stray cells or try slightly variating your resolution
133 std::cout << "Discrete Normal is not valid for interpolated pressure boundary!" << std::endl;
135 }
136 if (dynamics) {
137 dynamics->getParameters(_block).template set<descriptors::OMEGA>(omega);
138 }
139 setBoundary(_block, iX, iY, iZ, dynamics);
140 }
141 });
142}
143
144
145}
146
147#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.
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.
void exit(int exitcode)
Definition singleton.h:165
Top level namespace for all of OpenLB.
void setInterpolatedPressureBoundary(SuperLattice< T, DESCRIPTOR > &sLattice, T omega, SuperGeometry< T, 2 > &superGeometry, int material)
Initialising the setInterpolatedPressureBoundary 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