OpenLB 1.7
Loading...
Searching...
No Matches
Classes | Functions
olb::legacy Namespace Reference

Classes

struct  BasicDynamics
 
class  BGKdynamics
 
class  NoLatticeDynamics
 
class  OffDynamics
 Dynamics for offLattice boundary conditions OffDynamics are basically NoLatticeDynamics with the additional functionality to store given velocities exactly at boundary links. More...
 

Functions

template<typename T , typename DESCRIPTOR >
void defineUBouzidi (SuperLattice< T, DESCRIPTOR > &sLattice, SuperGeometry< T, 2 > &superGeometry, int material, AnalyticalF2D< T, T > &u, std::vector< int > bulkMaterials=std::vector< int >(1, 1))
 
template<typename T , typename DESCRIPTOR >
void defineUBouzidi (SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF2D< T > > &&indicator, FunctorPtr< SuperIndicatorF2D< T > > &&bulkIndicator, AnalyticalF2D< T, T > &u)
 
template<typename T , typename DESCRIPTOR >
void defineUBouzidi (BlockLattice< T, DESCRIPTOR > &block, BlockIndicatorF2D< T > &indicator, BlockIndicatorF2D< T > &bulkIndicator, AnalyticalF2D< T, T > &u)
 
template<typename T , typename DESCRIPTOR >
void defineUBouzidi (BlockLattice< T, DESCRIPTOR > &block, int iX, int iY, int iPop, const T u[DESCRIPTOR::d])
 
template<typename T , typename DESCRIPTOR >
bool getBoundaryIntersection (BlockLattice< T, DESCRIPTOR > &block, int iX, int iY, int iPop, T point[DESCRIPTOR::d])
 
template<typename T , typename DESCRIPTOR >
void setBoundaryIntersection (BlockLattice< T, DESCRIPTOR > &block, int iX, int iY, int iPop, T distance)
 
template<typename T , typename DESCRIPTOR >
void defineUBouzidi (SuperLattice< T, DESCRIPTOR > &sLattice, SuperGeometry< T, 3 > &superGeometry, int material, AnalyticalF3D< T, T > &u, std::vector< int > bulkMaterials=std::vector< int >(1, 1))
 
template<typename T , typename DESCRIPTOR >
void defineUBouzidi (SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF3D< T > > &&boundaryIndicator, AnalyticalF3D< T, T > &u, std::vector< int > bulkMaterials=std::vector< int >(1, 1))
 
template<typename T , typename DESCRIPTOR >
void defineUBouzidi (SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF3D< T > > &&boundaryIndicator, FunctorPtr< SuperIndicatorF3D< T > > &&bulkIndicator, AnalyticalF3D< T, T > &u)
 
template<typename T , typename DESCRIPTOR >
void defineUBouzidi (BlockLattice< T, DESCRIPTOR > &_block, BlockIndicatorF3D< T > &indicator, BlockIndicatorF3D< T > &bulkIndicator, AnalyticalF3D< T, T > &u)
 
template<typename T , typename DESCRIPTOR >
void defineUBouzidi (BlockLattice< T, DESCRIPTOR > &_block, int iX, int iY, int iZ, int iPop, const T u[DESCRIPTOR::d])
 
template<typename T , typename DESCRIPTOR >
bool getBoundaryIntersection (BlockLattice< T, DESCRIPTOR > &_block, int iX, int iY, int iZ, int iPop, T point[DESCRIPTOR::d])
 
template<typename T , typename DESCRIPTOR >
void setBoundaryIntersection (BlockLattice< T, DESCRIPTOR > &_block, int iX, int iY, int iZ, int iPop, T distance)
 
template<typename T , typename DESCRIPTOR , class MixinDynamics = BGKdynamics<T,DESCRIPTOR>>
void setBouzidiVelocityBoundary (SuperLattice< T, DESCRIPTOR > &sLattice, SuperGeometry< T, 2 > &superGeometry, int material, IndicatorF2D< T > &geometryIndicator, std::vector< int > bulkMaterials=std::vector< int >(1, 1))
 Set offDynamics with boundary links and post processors using indicators.
 
template<typename T , typename DESCRIPTOR , class MixinDynamics >
void setBouzidiVelocityBoundary (SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF2D< T > > &&boundaryIndicator, FunctorPtr< SuperIndicatorF2D< T > > &&bulkIndicator, IndicatorF2D< T > &geometryIndicator)
 Initialising the BouzidiVelocityBoundary on the superLattice domain.
 
template<typename T , typename DESCRIPTOR , class MixinDynamics >
void setBouzidiVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockIndicatorF2D< T > &boundaryIndicator, BlockIndicatorF2D< T > &bulkIndicator, IndicatorF2D< T > &geometryIndicator)
 
template<typename T , typename DESCRIPTOR , class MixinDynamics >
void setBouzidiVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 2 > &blockGeometryStructure, int iX, int iY, BlockIndicatorF2D< T > &bulkIndicator, IndicatorF2D< T > &geometryIndicator)
 
template<typename T , typename DESCRIPTOR , class MixinDynamics >
void setBouzidiVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 2 > &blockGeometryStructure, int x, int y, T distances[DESCRIPTOR::q])
 
template<typename T , typename DESCRIPTOR , class MixinDynamics >
void setBouzidiVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 2 > &blockGeometryStructure, int x, int y, int iPop, T dist)
 
template<typename T , typename DESCRIPTOR , class MixinDynamics >
void setOffDynamics (BlockLattice< T, DESCRIPTOR > &block, int x, int y, T location[DESCRIPTOR::d], T distances[DESCRIPTOR::q])
 
template<typename T , typename DESCRIPTOR >
void setBouzidiVelocityBoundary (SuperLattice< T, DESCRIPTOR > &sLattice, SuperGeometry< T, 3 > &superGeometry, int material, IndicatorF3D< T > &indicator, std::vector< int > bulkMaterials=std::vector< int >(1, 1))
 Set offDynamics with boundary links and post processors using indicators.
 
template<typename T , typename DESCRIPTOR >
void setBouzidiVelocityBoundary (SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF3D< T > > &&boundaryIndicator, FunctorPtr< SuperIndicatorF3D< T > > &&bulkIndicator, IndicatorF3D< T > &geometryIndicator)
 Initialising the BouzidiVelocityBoundary on the superLattice domain.
 
template<typename T , typename DESCRIPTOR >
void setBouzidiVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockIndicatorF3D< T > &boundaryIndicator, BlockIndicatorF3D< T > &bulkIndicator, IndicatorF3D< T > &geometryIndicator, T _epsFraction)
 
template<typename T , typename DESCRIPTOR >
void setBouzidiVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 3 > &blockGeometryStructure, int iX, int iY, int iZ, IndicatorF3D< T > &geometryIndicator, BlockIndicatorF3D< T > &bulkIndicator, T _epsFraction)
 
template<typename T , typename DESCRIPTOR >
void setBouzidiVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 3 > &blockGeometryStructure, int x, int y, int z, T distances[DESCRIPTOR::q])
 
template<typename T , typename DESCRIPTOR >
void setBouzidiVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 3 > &blockGeometryStructure, int x, int y, int z, int iPop, T dist)
 
template<typename T , typename DESCRIPTOR , class MixinDynamics = BGKdynamics<T,DESCRIPTOR>>
void setBouzidiZeroVelocityBoundary (SuperLattice< T, DESCRIPTOR > &sLattice, SuperGeometry< T, 2 > &superGeometry, int material, IndicatorF2D< T > &geometryIndicator, std::vector< int > bulkMaterials=std::vector< int >(1, 1))
 Set offDynamics with boundary links and post processors using indicators.
 
template<typename T , typename DESCRIPTOR , class MixinDynamics = BGKdynamics<T,DESCRIPTOR>>
void setBouzidiZeroVelocityBoundary (SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF2D< T > > &&boundaryIndicator, FunctorPtr< SuperIndicatorF2D< T > > &&bulkIndicator, IndicatorF2D< T > &geometryIndicator)
 Initialising the BouzidiZeroVelocityBoundary on the superLattice domain.
 
template<typename T , typename DESCRIPTOR , class MixinDynamics >
void setBouzidiZeroVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockIndicatorF2D< T > &boundaryIndicator, BlockIndicatorF2D< T > &bulkIndicator, IndicatorF2D< T > &geometryIndicator)
 
template<typename T , typename DESCRIPTOR , class MixinDynamics >
void setBouzidiZeroVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 2 > &blockGeometryStructure, int iX, int iY, BlockIndicatorF2D< T > &bulkIndicator, IndicatorF2D< T > &geometryIndicator)
 
template<typename T , typename DESCRIPTOR , class MixinDynamics >
void setBouzidiZeroVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 2 > &blockGeometryStructure, int x, int y, T distances[DESCRIPTOR::q])
 
template<typename T , typename DESCRIPTOR , class MixinDynamics >
void setBouzidiZeroVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 2 > &blockGeometryStructure, int x, int y, int iPop, T dist)
 
template<typename T , typename DESCRIPTOR >
void setBouzidiZeroVelocityBoundary (SuperLattice< T, DESCRIPTOR > &sLattice, SuperGeometry< T, 3 > &superGeometry, int material, IndicatorF3D< T > &geometryIndicator, std::vector< int > bulkMaterials=std::vector< int >(1, 1))
 Set offDynamics with boundary links and post processors using indicators.
 
template<typename T , typename DESCRIPTOR >
void setBouzidiZeroVelocityBoundary (SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF3D< T > > &&boundaryIndicator, IndicatorF3D< T > &geometryIndicator, std::vector< int > bulkMaterials)
 
template<typename T , typename DESCRIPTOR >
void setBouzidiZeroVelocityBoundary (SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF3D< T > > &&boundaryIndicator, FunctorPtr< SuperIndicatorF3D< T > > &&bulkIndicator, IndicatorF3D< T > &geometryIndicator)
 
template<typename T , typename DESCRIPTOR >
void setBouzidiZeroVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockIndicatorF3D< T > &boundaryIndicator, BlockIndicatorF3D< T > &bulkIndicator, IndicatorF3D< T > &geometryIndicator, T _epsFraction)
 
template<typename T , typename DESCRIPTOR >
void setBouzidiZeroVelocityBoundary1 (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 3 > &blockGeometryStructure, int iX, int iY, int iZ, IndicatorF3D< T > &geometryIndicator, BlockIndicatorF3D< T > &bulkIndicator, T _epsFraction)
 
template<typename T , typename DESCRIPTOR >
void setBouzidiZeroVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 3 > &blockGeometryStructure, int x, int y, int z, T distances[DESCRIPTOR::q])
 
template<typename T , typename DESCRIPTOR >
void setBouzidiZeroVelocityBoundary (BlockLattice< T, DESCRIPTOR > &block, BlockGeometry< T, 3 > &blockGeometryStructure, int x, int y, int z, int iPop, T dist)
 

Function Documentation

◆ defineUBouzidi() [1/9]

template<typename T , typename DESCRIPTOR >
void olb::legacy::defineUBouzidi ( BlockLattice< T, DESCRIPTOR > & _block,
BlockIndicatorF3D< T > & indicator,
BlockIndicatorF3D< T > & bulkIndicator,
AnalyticalF3D< T, T > & u )

Definition at line 67 of file defineU3D.hh.

68{
69 _block.forSpatialLocations([&](auto iX, auto iY, auto iZ) {
70 if (indicator(iX,iY,iZ)) {
71 for (int q = 1; q < DESCRIPTOR::q ; ++q) {
72 // Get direction
73 const int iXn = iX + descriptors::c<DESCRIPTOR>(q,0);
74 const int iYn = iY + descriptors::c<DESCRIPTOR>(q,1);
75 const int iZn = iZ + descriptors::c<DESCRIPTOR>(q,2);
76 if (_block.isInside({iXn,iYn,iZn}) && bulkIndicator(iXn,iYn,iZn)) {
77 T intersection[3] = { };
78 const int opp = descriptors::opposite<DESCRIPTOR>(q);
79 if (getBoundaryIntersection<T,DESCRIPTOR>(_block, iX, iY, iZ, opp, intersection)) {
80 T vel[3]= { };
81 u(vel, intersection);
82 defineUBouzidi<T,DESCRIPTOR>(_block, iX, iY, iZ, opp, vel);
83 }
84 }
85 }
86 }
87 });
88}
void forSpatialLocations(F f) const
bool isInside(LatticeR< D > latticeR) const
Return whether location is valid.

References olb::BlockStructureD< D >::forSpatialLocations(), and olb::BlockStructureD< D >::isInside().

+ Here is the call graph for this function:

◆ defineUBouzidi() [2/9]

template<typename T , typename DESCRIPTOR >
void olb::legacy::defineUBouzidi ( BlockLattice< T, DESCRIPTOR > & _block,
int iX,
int iY,
int iZ,
int iPop,
const T u[DESCRIPTOR::d] )

Definition at line 93 of file defineU3D.hh.

94{
95 OstreamManager clout(std::cout, "defineUBouzidi");
96 bool _output = false;
97 static_cast<legacy::OffDynamics<T,DESCRIPTOR>*>(_block.getDynamics(iX, iY, iZ))->defineU(iPop, u);
98 if (_output) {
99 clout << "defineUBouzidi(" << iX << ", " << iY << ", " << iZ << " )" << std::endl;
100 }
101
102}
virtual Dynamics< T, DESCRIPTOR > * getDynamics(DynamicsPromise< T, DESCRIPTOR > &&)=0
Return pointer to dynamics yielded by promise.
class for marking output with some text
Dynamics for offLattice boundary conditions OffDynamics are basically NoLatticeDynamics with the addi...
Definition dynamics.h:145

References olb::BlockLattice< T, DESCRIPTOR >::getDynamics().

+ Here is the call graph for this function:

◆ defineUBouzidi() [3/9]

template<typename T , typename DESCRIPTOR >
void olb::legacy::defineUBouzidi ( BlockLattice< T, DESCRIPTOR > & block,
BlockIndicatorF2D< T > & indicator,
BlockIndicatorF2D< T > & bulkIndicator,
AnalyticalF2D< T, T > & u )

Definition at line 60 of file defineU2D.hh.

61{
62 block.forSpatialLocations([&](auto iX, auto iY) {
63 if (indicator(iX,iY)) {
64 for (int iPop = 1; iPop < DESCRIPTOR::q ; ++iPop) {
65 int iXn = iX + descriptors::c<DESCRIPTOR >(iPop,0);
66 int iYn = iY + descriptors::c<DESCRIPTOR >(iPop,1);
67 if (block.isInside({iXn,iYn}) && bulkIndicator(iXn, iYn)) {
68 T intersection[] = { T(), T() }; // coord. of intersection
69 int opp = descriptors::opposite<DESCRIPTOR >(iPop);
70 if (getBoundaryIntersection<T,DESCRIPTOR>(block, iX, iY, opp, intersection) ) {
71 T vel[]= {T(),T()};
72 u(vel,intersection);
73 defineUBouzidi<T,DESCRIPTOR>(iX, iY, opp, vel);
74 }
75 }
76 }
77 }
78 });
79}

References olb::BlockStructureD< D >::forSpatialLocations(), and olb::BlockStructureD< D >::isInside().

+ Here is the call graph for this function:

◆ defineUBouzidi() [4/9]

template<typename T , typename DESCRIPTOR >
void olb::legacy::defineUBouzidi ( BlockLattice< T, DESCRIPTOR > & block,
int iX,
int iY,
int iPop,
const T u[DESCRIPTOR::d] )

Definition at line 83 of file defineU2D.hh.

84{
85 bool _output = false;
86 OstreamManager clout(std::cout, "defineUBouzidi");
87 block.getDynamics(iX, iY)->defineU(iPop, u);
88 if (_output) {
89 clout << "defineUBouzidi(" << iX << ", " << iY << " )" << std::endl;
90 }
91}

References olb::BlockLattice< T, DESCRIPTOR >::getDynamics().

+ Here is the call graph for this function:

◆ defineUBouzidi() [5/9]

template<typename T , typename DESCRIPTOR >
void olb::legacy::defineUBouzidi ( SuperLattice< T, DESCRIPTOR > & sLattice,
FunctorPtr< SuperIndicatorF2D< T > > && indicator,
FunctorPtr< SuperIndicatorF2D< T > > && bulkIndicator,
AnalyticalF2D< T, T > & u )

Definition at line 45 of file defineU2D.hh.

48{
49
50 for (int iCloc = 0; iCloc < sLattice.getLoadBalancer().size(); ++iCloc) {
51 defineUBouzidi<T,DESCRIPTOR>(sLattice.getBlockIndicator(iCloc), indicator->getBlockIndicatorF(iCloc),
52 bulkIndicator->getBlockIndicatorF(iCloc),
53 u);
54 }
55}
LoadBalancer< T > & getLoadBalancer()
Read and write access to the load balancer.

References olb::SuperStructure< T, D >::getLoadBalancer().

+ Here is the call graph for this function:

◆ defineUBouzidi() [6/9]

template<typename T , typename DESCRIPTOR >
void olb::legacy::defineUBouzidi ( SuperLattice< T, DESCRIPTOR > & sLattice,
FunctorPtr< SuperIndicatorF3D< T > > && boundaryIndicator,
AnalyticalF3D< T, T > & u,
std::vector< int > bulkMaterials = std::vector<int>(1,1) )

Definition at line 43 of file defineU3D.hh.

45{
46 defineUBouzidi<T,DESCRIPTOR>(sLattice, std::forward<decltype(boundaryIndicator)>(boundaryIndicator),
47 boundaryIndicator->getSuperGeometry().getMaterialIndicator(std::move(bulkMaterials)),
48 u);
49}

◆ defineUBouzidi() [7/9]

template<typename T , typename DESCRIPTOR >
void olb::legacy::defineUBouzidi ( SuperLattice< T, DESCRIPTOR > & sLattice,
FunctorPtr< SuperIndicatorF3D< T > > && boundaryIndicator,
FunctorPtr< SuperIndicatorF3D< T > > && bulkIndicator,
AnalyticalF3D< T, T > & u )

Definition at line 52 of file defineU3D.hh.

54{
55
56 for (int iCloc = 0; iCloc < sLattice.getLoadBalancer().size(); ++iCloc) {
57 defineUBouzidi<T,DESCRIPTOR>(sLattice.getBlock(iCloc),
58 boundaryIndicator->getBlockIndicatorF(iCloc),
59 bulkIndicator->getBlockIndicatorF(iCloc),u);
60 }
61
62}
BlockLattice< T, DESCRIPTOR > & getBlock(int locC)
Return BlockLattice with local index locC.

References olb::SuperLattice< T, DESCRIPTOR >::getBlock(), and olb::SuperStructure< T, D >::getLoadBalancer().

+ Here is the call graph for this function:

◆ defineUBouzidi() [8/9]

template<typename T , typename DESCRIPTOR >
void olb::legacy::defineUBouzidi ( SuperLattice< T, DESCRIPTOR > & sLattice,
SuperGeometry< T, 2 > & superGeometry,
int material,
AnalyticalF2D< T, T > & u,
std::vector< int > bulkMaterials = std::vector<int>(1,1) )

Definition at line 35 of file defineU2D.hh.

38{
39 defineUBouzidi<T,DESCRIPTOR>(sLattice, superGeometry.getMaterialIndicator(material),
40 superGeometry.getMaterialIndicator(std::move(bulkMaterials)),
41 u);
42}
std::unique_ptr< SuperIndicatorF< T, D > > getMaterialIndicator(std::vector< int > &&materials)
Returns a material indicator using the given vector of materials.

References olb::SuperGeometry< T, D >::getMaterialIndicator().

+ Here is the call graph for this function:

◆ defineUBouzidi() [9/9]

template<typename T , typename DESCRIPTOR >
void olb::legacy::defineUBouzidi ( SuperLattice< T, DESCRIPTOR > & sLattice,
SuperGeometry< T, 3 > & superGeometry,
int material,
AnalyticalF3D< T, T > & u,
std::vector< int > bulkMaterials = std::vector<int>(1,1) )

Definition at line 35 of file defineU3D.hh.

37{
38 defineUBouzidi<T,DESCRIPTOR>(sLattice, superGeometry.getMaterialIndicator(material), u, bulkMaterials);
39}

References olb::SuperGeometry< T, D >::getMaterialIndicator().

+ Here is the call graph for this function:

◆ getBoundaryIntersection() [1/2]

template<typename T , typename DESCRIPTOR >
bool olb::legacy::getBoundaryIntersection ( BlockLattice< T, DESCRIPTOR > & _block,
int iX,
int iY,
int iZ,
int iPop,
T point[DESCRIPTOR::d] )

Definition at line 105 of file defineU3D.hh.

106{
107 return static_cast<legacy::OffDynamics<T,DESCRIPTOR>*>(_block.getDynamics(iX, iY, iZ))->getBoundaryIntersection(iPop, point);
108}
bool getBoundaryIntersection(BlockLattice< T, DESCRIPTOR > &block, int iX, int iY, int iPop, T point[DESCRIPTOR::d])
Definition defineU2D.hh:96

References getBoundaryIntersection(), and olb::BlockLattice< T, DESCRIPTOR >::getDynamics().

+ Here is the call graph for this function:

◆ getBoundaryIntersection() [2/2]

template<typename T , typename DESCRIPTOR >
bool olb::legacy::getBoundaryIntersection ( BlockLattice< T, DESCRIPTOR > & block,
int iX,
int iY,
int iPop,
T point[DESCRIPTOR::d] )

Definition at line 96 of file defineU2D.hh.

97{
98 return static_cast<legacy::OffDynamics<T,DESCRIPTOR>*>(block.getDynamics(iX,iY))->getBoundaryIntersection(iPop, point);
99}

References getBoundaryIntersection(), and olb::BlockLattice< T, DESCRIPTOR >::getDynamics().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ setBoundaryIntersection() [1/2]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBoundaryIntersection ( BlockLattice< T, DESCRIPTOR > & _block,
int iX,
int iY,
int iZ,
int iPop,
T distance )

Definition at line 111 of file defineU3D.hh.

112{
113 bool _output = false;
114 OstreamManager clout(std::cout, "setBoundaryIntersection");
115 _block.getDynamics(iX, iY, iZ)->setBoundaryIntersection(iPop, distance);
116 if (_output) {
117 clout << "setBoundaryIntersection(" << iX << ", " << iY << ", " << iZ << " )" << std::endl;
118 }
119}

References olb::BlockLattice< T, DESCRIPTOR >::getDynamics().

+ Here is the call graph for this function:

◆ setBoundaryIntersection() [2/2]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBoundaryIntersection ( BlockLattice< T, DESCRIPTOR > & block,
int iX,
int iY,
int iPop,
T distance )

Definition at line 102 of file defineU2D.hh.

103{
104 bool _output = false;
105 OstreamManager clout(std::cout, "setBoundaryIntersection");
106 static_cast<legacy::OffDynamics<T,DESCRIPTOR>*>(block.getDynamics(iX,iY))->setBoundaryIntersection(iPop, distance);
107 if (_output) {
108 clout << "setBoundaryIntersection(" << iX << ", " << iY << " )" << std::endl;
109 }
110}
void setBoundaryIntersection(BlockLattice< T, DESCRIPTOR > &block, int iX, int iY, int iPop, T distance)
Definition defineU2D.hh:102

References olb::BlockLattice< T, DESCRIPTOR >::getDynamics(), and setBoundaryIntersection().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ setBouzidiVelocityBoundary() [1/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics >
void olb::legacy::setBouzidiVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 2 > & blockGeometryStructure,
int iX,
int iY,
BlockIndicatorF2D< T > & bulkIndicator,
IndicatorF2D< T > & geometryIndicator )

Definition at line 101 of file setBouzidiVelocityBoundary2D.hh.

102{
103 T _epsFraction = 0.0001;
104 OstreamManager clout(std::cout, "setBouzidiVelocityBoundary");
105 T distances[DESCRIPTOR::q];
106 for (int iPop = 1; iPop < DESCRIPTOR::q ; ++iPop) {
107 distances[iPop] = -1;
108 }
109
110 for (int iPop = 1; iPop < DESCRIPTOR::q ; ++iPop) {
111 const int iXn = iX + descriptors::c<DESCRIPTOR >(iPop,0);
112 const int iYn = iY + descriptors::c<DESCRIPTOR >(iPop,1);
113 if (blockGeometryStructure.isInside(iXn,iYn) && bulkIndicator(iXn,iYn)) {
114 T dist = -1;
115 T physR[2];
116 blockGeometryStructure.getPhysR(physR,{iXn,iYn});
117 T voxelSize=blockGeometryStructure.getDeltaR();
118 Vector<T,2> physC(physR);
119
120 Vector<T,2> direction(-voxelSize*descriptors::c<DESCRIPTOR >(iPop,0),-voxelSize*descriptors::c<DESCRIPTOR >(iPop,1));
121 T cPhysNorm = voxelSize*util::sqrt(descriptors::c<DESCRIPTOR >(iPop,0)*descriptors::c<DESCRIPTOR >(iPop,0)+descriptors::c<DESCRIPTOR >(iPop,1)*descriptors::c<DESCRIPTOR >(iPop,1));
122
123 if (!geometryIndicator.distance(dist,physC,direction,blockGeometryStructure.getIcGlob() ) ) {
124 T epsX = voxelSize*descriptors::c<DESCRIPTOR >(iPop,0)*_epsFraction;
125 T epsY = voxelSize*descriptors::c<DESCRIPTOR >(iPop,1)*_epsFraction;
126
127 Vector<T,2> physC2(physC);
128 physC2[0] += epsX;
129 physC2[1] += epsY;
130 Vector<T,2> direction2(direction);
131 direction2[0] -= 2.*epsX;
132 direction2[1] -= 2.*epsY;
133
134 if ( !geometryIndicator.distance(dist,physC2,direction2,blockGeometryStructure.getIcGlob())) {
135 clout << "ERROR: no boundary found at (" << iXn << "," << iYn <<") ~ ("
136 << physR[0] << "," << physR[1] << "), "
137 << "in direction " << descriptors::opposite<DESCRIPTOR >(iPop)
138 << std::endl;
139 }
140 T distNew = (dist - util::sqrt(epsX*epsX+epsY*epsY))/cPhysNorm;
141 if (distNew < 0.5) {
142 dist = 0;
143 }
144 else {
145 dist = 0.5 * cPhysNorm;
146 clout << "WARNING: distance at (" << iXn << "," << iYn <<") ~ ("
147 << physR[0] << "," << physR[1] <<"), "
148 << "in direction " << descriptors::opposite<DESCRIPTOR >(iPop) << ": "
149 << distNew
150 << " rounded to "
151 << dist/cPhysNorm
152 << std::endl;
153 }
154 }
155 distances[descriptors::opposite<DESCRIPTOR >(iPop)] = dist/cPhysNorm;
156 } // bulk indicator
157 } // iPop loop
158 setBouzidiVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(block, blockGeometryStructure, iX, iY, distances);
159
160}
Vector< T, D > getPhysR(LatticeR< D > latticeR)
int getIcGlob() const
Read only access to the global iC number which is given !=-1 if the block geometries are part of a su...
T getDeltaR() const
Read only access to the voxel size given in SI units (meter)
virtual bool distance(S &distance, const Vector< S, 2 > &origin, S precision, const Vector< S, 2 > &direction)
returns false or true and pos. distance if there was one found for an given origin and direction

References olb::IndicatorF2D< S >::distance(), olb::BlockGeometry< T, D >::getDeltaR(), olb::BlockGeometry< T, D >::getIcGlob(), olb::BlockGeometry< T, D >::getPhysR(), olb::BlockStructureD< D >::isInside(), and olb::util::sqrt().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [2/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics >
void olb::legacy::setBouzidiVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 2 > & blockGeometryStructure,
int x,
int y,
int iPop,
T dist )

Using Bouzidi BC OnePoint. This corresponds to Bounce Back

Using Bouzidi BC TwoPoint. This corresponds to Linear Interpolation

Using Bouzidi BC OnePoint. This corresponds to Bounce Back

Using Bouzidi BC TwoPoint. This corresponds to Linear Interpolation

Definition at line 190 of file setBouzidiVelocityBoundary2D.hh.

191{
193 if (blockGeometryStructure.getMaterial({x-descriptors::c<DESCRIPTOR >(iPop,0), y-descriptors::c<DESCRIPTOR >(iPop,1)}) != 1) {
195 postProcessor = new VelocityBounceBackPostProcessorGenerator2D<T,DESCRIPTOR>(x, y, iPop, dist);
196
197 }
198 else {
200 postProcessor =new VelocityBouzidiLinearPostProcessorGenerator2D<T,DESCRIPTOR>(x, y, iPop, dist);
201 }
202
203 if (postProcessor && !block.isPadding({x,y})) {
204 block.addPostProcessor(*postProcessor);
205 }
206}
int getMaterial(LatticeR< D > latticeR) const
returns the (iX,iY) entry in the 2D scalar field
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.
bool isPadding(LatticeR< D > latticeR) const
Return whether location is valid.

References olb::BlockLattice< T, DESCRIPTOR >::addPostProcessor(), olb::BlockGeometry< T, D >::getMaterial(), and olb::BlockStructureD< D >::isPadding().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [3/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics >
void olb::legacy::setBouzidiVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 2 > & blockGeometryStructure,
int x,
int y,
T distances[DESCRIPTOR::q] )

Definition at line 163 of file setBouzidiVelocityBoundary2D.hh.

164{
165 typedef DESCRIPTOR L;
166 T location[DESCRIPTOR::d];
167 blockGeometryStructure.getPhysR(location, {x,y});
168
169 T distancesCopy[L::q];
170 T spacing = blockGeometryStructure.getDeltaR();
171 for (int iPop = 1; iPop < L::q ; ++iPop) {
172 distancesCopy[iPop] = spacing*(1.-distances[iPop]);
173 if ( !util::nearZero(distances[iPop]+1) ) {
174 distancesCopy[iPop] = -1;
175 }
176 }
177 //can be universally used for Bouzidi and bounceBackVelocityBoundary
178 setOffDynamics<T,DESCRIPTOR,MixinDynamics>(block, x, y, location, distancesCopy);
179
180 for (int iPop = 1; iPop < L::q ; ++iPop) {
181 if ( !util::nearZero(distances[iPop]+1) ) {
182 setBouzidiVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(block, blockGeometryStructure, x-descriptors::c<L>(iPop,0), y-descriptors::c<L>(iPop,1), iPop, distances[iPop]);
183 }
184 }
185
186}

References olb::BlockGeometry< T, D >::getDeltaR(), olb::BlockGeometry< T, D >::getPhysR(), and olb::util::nearZero().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [4/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 3 > & blockGeometryStructure,
int iX,
int iY,
int iZ,
IndicatorF3D< T > & geometryIndicator,
BlockIndicatorF3D< T > & bulkIndicator,
T _epsFraction )

Definition at line 108 of file setBouzidiVelocityBoundary3D.hh.

110{
111 OstreamManager clout(std::cout, "setBouzidiBoundary");
112 T distances[DESCRIPTOR::q];
113 for (int iPop = 1; iPop < DESCRIPTOR::q ; ++iPop) {
114 distances[iPop] = -1;
115 }
116
117 for (int iPop = 1; iPop < DESCRIPTOR::q ; ++iPop) {
118 const Vector<int,3> c = descriptors::c<DESCRIPTOR>(iPop);
119 const int iXn = iX + c[0];
120 const int iYn = iY + c[1];
121 const int iZn = iZ + c[2];
122 if (blockGeometryStructure.isInside(iXn,iYn,iZn) && bulkIndicator(iXn,iYn,iZn)) {
123 T dist = -1;
124 T physR[3];
125 blockGeometryStructure.getPhysR(physR,{iXn,iYn,iZn});
126 T voxelSize=blockGeometryStructure.getDeltaR();
127
128 Vector<T,3> physC(physR);
129 Vector<T,3> direction(-voxelSize*c[0],-voxelSize*c[1],-voxelSize*c[2]);
130 T cPhysNorm = voxelSize*util::sqrt(c[0]*c[0]+c[1]*c[1]+c[2]*c[2]);
131
132 if (!geometryIndicator.distance(dist,physC,direction,blockGeometryStructure.getIcGlob() ) ) {
133 T epsX = voxelSize*c[0]*_epsFraction;
134 T epsY = voxelSize*c[1]*_epsFraction;
135 T epsZ = voxelSize*c[2]*_epsFraction;
136
137 Vector<T,3> physC2(physC);
138 physC2[0] += epsX;
139 physC2[1] += epsY;
140 physC2[2] += epsZ;
141 Vector<T,3> direction2(direction);
142 direction2[0] -= 2.*epsX;
143 direction2[1] -= 2.*epsY;
144 direction2[2] -= 2.*epsZ;
145
146 if ( !geometryIndicator.distance(dist,physC2,direction2,blockGeometryStructure.getIcGlob())) {
147 clout << "ERROR: no boundary found at (" << iXn << "," << iYn << "," << iZn <<") ~ ("
148 << physR[0] << "," << physR[1] << "," << physR[2] <<"), "
149 << "in direction " << descriptors::opposite<DESCRIPTOR >(iPop)
150 << std::endl;
151
152 }
153 T distNew = (dist - util::sqrt(epsX*epsX+epsY*epsY+epsZ*epsZ))/cPhysNorm;
154 if (distNew < 0.5) {
155 dist = 0;
156 }
157 else {
158 dist = 0.5 * cPhysNorm;
159 clout << "WARNING: distance at (" << iXn << "," << iYn << "," << iZn <<") ~ ("
160 << physR[0] << "," << physR[1] << "," << physR[2] <<"), "
161 << "in direction " << descriptors::opposite<DESCRIPTOR >(iPop) << ": "
162 << distNew
163 << " rounded to "
164 << dist/cPhysNorm
165 << std::endl;
166
167 }
168 }
169 distances[descriptors::opposite<DESCRIPTOR >(iPop)] = dist/cPhysNorm;
170 } // bulk indicator if
171 } // iPop loop
172 setBouzidiVelocityBoundary<T,DESCRIPTOR>(block, blockGeometryStructure, iX, iY, iZ, distances);
173}
virtual bool distance(S &distance, const Vector< S, 3 > &origin, S precision, const Vector< S, 3 > &direction)
Plain old scalar vector.
Definition vector.h:47

References olb::IndicatorF3D< S >::distance(), olb::BlockGeometry< T, D >::getDeltaR(), olb::BlockGeometry< T, D >::getIcGlob(), olb::BlockGeometry< T, D >::getPhysR(), olb::BlockStructureD< D >::isInside(), and olb::util::sqrt().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [5/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 3 > & blockGeometryStructure,
int x,
int y,
int z,
int iPop,
T dist )

Definition at line 202 of file setBouzidiVelocityBoundary3D.hh.

203{
204 const Vector<int,3> c = descriptors::c<DESCRIPTOR>(iPop);
205 //setOnePointVelocityBoundary. This corresponds to Bounce Back
206 auto postProcessor = std::unique_ptr<PostProcessorGenerator3D<T, DESCRIPTOR>> { nullptr };
207 if (blockGeometryStructure.getMaterial({x-c[0], y-c[1], z-c[2]}) != 1) {
208 postProcessor.reset(new VelocityBounceBackPostProcessorGenerator3D <T, DESCRIPTOR>(x, y, z, iPop, dist));
209 }
210 //setTwoPointVelocityBoundary.This corresponds to Linear Interpolation
211 else {
212 postProcessor.reset(new VelocityBouzidiLinearPostProcessorGenerator3D<T, DESCRIPTOR>(x, y, z, iPop, dist));
213 }
214 if (postProcessor && !block.isPadding({x,y,z})) {
215 block.addPostProcessor(*postProcessor);
216 }
217}

References olb::BlockLattice< T, DESCRIPTOR >::addPostProcessor(), olb::BlockGeometry< T, D >::getMaterial(), and olb::BlockStructureD< D >::isPadding().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [6/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 3 > & blockGeometryStructure,
int x,
int y,
int z,
T distances[DESCRIPTOR::q] )

Definition at line 177 of file setBouzidiVelocityBoundary3D.hh.

178{
179 T location[DESCRIPTOR::d];
180 blockGeometryStructure.getPhysR(location, {x,y,z});
181 T distancesCopy[DESCRIPTOR::q];
182 T spacing = blockGeometryStructure.getDeltaR();
183 for (int iPop = 1; iPop < DESCRIPTOR::q ; ++iPop) {
184 distancesCopy[iPop] = spacing*(1.-distances[iPop]);
185 if ( util::nearZero(distances[iPop]+1) ) {
186 distancesCopy[iPop] = -1;
187 }
188 }
189 Dynamics<T,DESCRIPTOR>* dynamics = new legacy::OffDynamics<T, DESCRIPTOR>(location, distancesCopy);
190 block.defineDynamics({x,y,z}, dynamics);
191
192 for (int iPop = 1; iPop < DESCRIPTOR::q ; ++iPop) {
193 if (!util::nearZero(distances[iPop]+1)) {
194 const Vector<int,3> c = descriptors::c<DESCRIPTOR>(iPop);
195
196 setBouzidiVelocityBoundary<T,DESCRIPTOR>(block, blockGeometryStructure, x-c[0], y-c[1], z-c[2], iPop, distances[iPop]);
197 }
198 }
199}
void defineDynamics(LatticeR< DESCRIPTOR::d > latticeR, DynamicsPromise< T, DESCRIPTOR > &&promise)
Assign promised DYNAMICS to latticeR.
platform_constant int c[Q][D]

References olb::BlockLattice< T, DESCRIPTOR >::defineDynamics(), olb::BlockGeometry< T, D >::getDeltaR(), olb::BlockGeometry< T, D >::getPhysR(), and olb::util::nearZero().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [7/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics >
void olb::legacy::setBouzidiVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockIndicatorF2D< T > & boundaryIndicator,
BlockIndicatorF2D< T > & bulkIndicator,
IndicatorF2D< T > & geometryIndicator )

Definition at line 87 of file setBouzidiVelocityBoundary2D.hh.

88{
89 block.forSpatialLocations([&](auto iX, auto iY) {
90 if (boundaryIndicator(iX,iY)) {
91 setBouzidiVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(block,
92 bulkIndicator.getBlockGeometry(),
93 iX, iY,
94 bulkIndicator,
95 geometryIndicator);
96 }
97 });
98}
BlockGeometry< T, 2 > & getBlockGeometry()
Get underlying block geometry structure.

References olb::BlockStructureD< D >::forSpatialLocations(), and olb::BlockIndicatorF2D< T >::getBlockGeometry().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [8/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockIndicatorF3D< T > & boundaryIndicator,
BlockIndicatorF3D< T > & bulkIndicator,
IndicatorF3D< T > & geometryIndicator,
T _epsFraction )

Definition at line 91 of file setBouzidiVelocityBoundary3D.hh.

93{
94 block.forSpatialLocations([&](auto iX, auto iY, auto iZ) {
95 if (boundaryIndicator(iX,iY,iZ)) {
96 setBouzidiVelocityBoundary<T,DESCRIPTOR>(block,
97 bulkIndicator.getBlockGeometry(),
98 iX, iY, iZ,
99 geometryIndicator,
100 bulkIndicator,
101 _epsFraction);
102 }
103 });
104}
BlockGeometry< T, 3 > & getBlockGeometry()
Get underlying block geometry structure.

References olb::BlockStructureD< D >::forSpatialLocations(), and olb::BlockIndicatorF3D< T >::getBlockGeometry().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [9/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics >
void olb::legacy::setBouzidiVelocityBoundary ( SuperLattice< T, DESCRIPTOR > & sLattice,
FunctorPtr< SuperIndicatorF2D< T > > && boundaryIndicator,
FunctorPtr< SuperIndicatorF2D< T > > && bulkIndicator,
IndicatorF2D< T > & geometryIndicator )

Initialising the BouzidiVelocityBoundary on the superLattice domain.

Adds needed Cells to the Communicator _commBC in SuperLattice

Adds needed Cells to the Communicator _commBC in SuperLattice

Definition at line 52 of file setBouzidiVelocityBoundary2D.hh.

55{
56 T _epsFraction = 0.0001;
57 /* local boundaries: _overlap = 0;
58 * interp boundaries: _overlap = 1;
59 * bouzidi boundaries: _overlap = 1;
60 * extField boundaries: _overlap = 1;
61 * advectionDiffusion boundaries: _overlap = 1;
62 */
63 int _overlap = 1;
64
65 OstreamManager clout(std::cout, "setBouzidiVelocityBoundary");
66 clout << "epsFraction=" << _epsFraction << std::endl;
67 clout.setMultiOutput(true);
68 for (int iCloc = 0; iCloc < sLattice.getLoadBalancer().size(); ++iCloc) {
69 clout << "Cuboid globiC " << sLattice.getLoadBalancer().glob(iCloc)
70 << " starts to read distances for Velocity Boundary..." << std::endl;
71 setBouzidiVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice.getBlock(iCloc),boundaryIndicator->getBlockIndicatorF(iCloc),
72 bulkIndicator->getBlockIndicatorF(iCloc),
73 geometryIndicator);
74 clout << "Cuboid globiC " << sLattice.getLoadBalancer().glob(iCloc)
75 << " finished reading distances for Velocity Boundary." << std::endl;
76 }
77 clout.setMultiOutput(false);
79 addPoints2CommBC<T,DESCRIPTOR>(sLattice, std::forward<decltype(boundaryIndicator)>(boundaryIndicator), _overlap);
80}

References olb::SuperLattice< T, DESCRIPTOR >::getBlock(), olb::SuperStructure< T, D >::getLoadBalancer(), and olb::OstreamManager::setMultiOutput().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [10/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiVelocityBoundary ( SuperLattice< T, DESCRIPTOR > & sLattice,
FunctorPtr< SuperIndicatorF3D< T > > && boundaryIndicator,
FunctorPtr< SuperIndicatorF3D< T > > && bulkIndicator,
IndicatorF3D< T > & geometryIndicator )

Initialising the BouzidiVelocityBoundary on the superLattice domain.

Adds needed Cells to the Communicator _commBC in SuperLattice

Adds needed Cells to the Communicator _commBC in SuperLattice

Definition at line 52 of file setBouzidiVelocityBoundary3D.hh.

55{
56 OstreamManager clout(std::cout, "setBouzidiBoundary");
57 T _epsFraction = 0.0001;
58 int overlap = 1;
59 bool _output = false;
60 if (_output) {
61 clout << "epsFraction=" << _epsFraction << std::endl;
62 clout.setMultiOutput(true);
63 }
64 for (int iCloc = 0; iCloc < sLattice.getLoadBalancer().size(); ++iCloc) {
65 if (_output) {
66 clout << "Cuboid globiC " << sLattice.getLoadBalancer().glob(iCloc)
67 << " starts to read distances for Velocity Boundary..." << std::endl;
68 }
69 //this is a guess to solve the private member block issue
70 setBouzidiVelocityBoundary<T,DESCRIPTOR>(
71 sLattice.getBlock(iCloc), boundaryIndicator->getBlockIndicatorF(iCloc),
72 bulkIndicator->getBlockIndicatorF(iCloc),
73 geometryIndicator, _epsFraction);
74 if (_output) {
75 clout << "Cuboid globiC " << sLattice.getLoadBalancer().glob(iCloc)
76 << " finished reading distances for Velocity Boundary." << std::endl;
77 }
78 }
79 if (_output) {
80 clout.setMultiOutput(false);
81 }
83 addPoints2CommBC<T,DESCRIPTOR>(sLattice,std::forward<decltype(boundaryIndicator)>(boundaryIndicator), overlap);
84}

References olb::SuperLattice< T, DESCRIPTOR >::getBlock(), olb::SuperStructure< T, D >::getLoadBalancer(), and olb::OstreamManager::setMultiOutput().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [11/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics = BGKdynamics<T,DESCRIPTOR>>
void olb::legacy::setBouzidiVelocityBoundary ( SuperLattice< T, DESCRIPTOR > & sLattice,
SuperGeometry< T, 2 > & superGeometry,
int material,
IndicatorF2D< T > & geometryIndicator,
std::vector< int > bulkMaterials = std::vector<int>(1,1) )

Set offDynamics with boundary links and post processors using indicators.

Initialising the BouzidiVelocityBoundary on the superLattice domain.

Add offDynamics with initialisation of boundary links and the corresponding post processors Note: Uses information of the second neighbours of the cell (x,y) Add post processors. Ensure that offDynamics are defined!

Parameters
boundaryIndicatorIndicator describing boundary cells
bulkIndicatorIndicator describing bulk cells
geometryIndicatorIndicator describing the geometry to be bounded Initialising the BouzidiVelocityBoundary on the superLattice domain

Definition at line 41 of file setBouzidiVelocityBoundary2D.hh.

43{
44
45 setBouzidiVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice, superGeometry.getMaterialIndicator(material),
46 superGeometry.getMaterialIndicator(std::move(bulkMaterials)),
47 geometryIndicator);
48}

References olb::SuperGeometry< T, D >::getMaterialIndicator().

+ Here is the call graph for this function:

◆ setBouzidiVelocityBoundary() [12/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiVelocityBoundary ( SuperLattice< T, DESCRIPTOR > & sLattice,
SuperGeometry< T, 3 > & superGeometry,
int material,
IndicatorF3D< T > & indicator,
std::vector< int > bulkMaterials = std::vector<int>(1,1) )

Set offDynamics with boundary links and post processors using indicators.

Initialising the BouzidiVelocityBoundary on the superLattice domain.

Add offDynamics with initialisation of boundary links and the corresponding post processors Note: Uses information of the second neighbours of the cell (x,y,z) Add post processors. Ensure that offDynamics are defined!

Parameters
boundaryIndicatorIndicator describing boundary cells
bulkIndicatorIndicator describing bulk cells
geometryIndicatorIndicator describing the geometry to be bounded Initialising the BouzidiVelocityBoundary on the superLattice domain

Definition at line 40 of file setBouzidiVelocityBoundary3D.hh.

43{
44
45 setBouzidiVelocityBoundary<T,DESCRIPTOR>(sLattice, superGeometry.getMaterialIndicator(material),
46 superGeometry.getMaterialIndicator(std::move(bulkMaterials)),
47 geometryIndicator);
48}

References olb::SuperGeometry< T, D >::getMaterialIndicator().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [1/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics >
void olb::legacy::setBouzidiZeroVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 2 > & blockGeometryStructure,
int iX,
int iY,
BlockIndicatorF2D< T > & bulkIndicator,
IndicatorF2D< T > & geometryIndicator )

Definition at line 99 of file setBouzidiZeroVelocityBoundary2D.hh.

100{
101 OstreamManager clout(std::cout, "setBouzidiZeroVelocityBoundary");
102 T _epsFraction = 0.0001;
103 T distances[DESCRIPTOR::q];
104 std::fill(std::begin(distances), std::end(distances), -1);
105
106 for (int iPop = 1; iPop < DESCRIPTOR::q ; ++iPop) {
107 const Vector<int,2> c = descriptors::c<DESCRIPTOR>(iPop);
108 const int iXn = iX + c[0];
109 const int iYn = iY + c[1];
110 if (blockGeometryStructure.isInside(iXn,iYn)) {
111 if (bulkIndicator(iXn,iYn)) {
112 T dist = -1;
113 T physR[2];
114 blockGeometryStructure.getPhysR(physR,{iXn,iYn});
115 T voxelSize=blockGeometryStructure.getDeltaR();
116 Vector<T,2> physC(physR);
117
118 Vector<T,2> direction(-voxelSize*descriptors::c<DESCRIPTOR >(iPop,0),-voxelSize*descriptors::c<DESCRIPTOR >(iPop,1));
119 T cPhysNorm = voxelSize*util::sqrt(descriptors::c<DESCRIPTOR >(iPop,0)*descriptors::c<DESCRIPTOR >(iPop,0)+descriptors::c<DESCRIPTOR >(iPop,1)*descriptors::c<DESCRIPTOR >(iPop,1));
120
121 if (!geometryIndicator.distance(dist,physC,direction,blockGeometryStructure.getIcGlob() ) ) {
122 T epsX = voxelSize*descriptors::c<DESCRIPTOR >(iPop,0)*_epsFraction;
123 T epsY = voxelSize*descriptors::c<DESCRIPTOR >(iPop,1)*_epsFraction;
124
125 Vector<T,2> physC2(physC);
126 physC2[0] += epsX;
127 physC2[1] += epsY;
128 Vector<T,2> direction2(direction);
129 direction2[0] -= 2.*epsX;
130 direction2[1] -= 2.*epsY;
131
132 if ( !geometryIndicator.distance(dist,physC2,direction2,blockGeometryStructure.getIcGlob())) {
133 clout << "ERROR: no boundary found at (" << iXn << "," << iYn <<") ~ ("
134 << physR[0] << "," << physR[1] << "), "
135 << "in direction " << descriptors::opposite<DESCRIPTOR >(iPop)
136 << std::endl;
137 }
138 T distNew = (dist - util::sqrt(epsX*epsX+epsY*epsY))/cPhysNorm;
139 if (distNew < 0.5) {
140 dist = 0;
141 }
142 else {
143 dist = 0.5 * cPhysNorm;
144 clout << "WARNING: distance at (" << iXn << "," << iYn <<") ~ ("
145 << physR[0] << "," << physR[1] <<"), "
146 << "in direction " << descriptors::opposite<DESCRIPTOR >(iPop) << ": "
147 << distNew
148 << " rounded to "
149 << dist/cPhysNorm
150 << std::endl;
151 }
152 }
153 distances[descriptors::opposite<DESCRIPTOR >(iPop)] = dist/cPhysNorm;
154 }
155 else {
156 if (blockGeometryStructure.getMaterial({iXn,iYn}) != 0) {
157 auto postProcessor = std::unique_ptr<ZeroVelocityBounceBackPostProcessorGenerator2D<T,DESCRIPTOR>>{new ZeroVelocityBounceBackPostProcessorGenerator2D<T,DESCRIPTOR>(iXn, iYn, descriptors::opposite<DESCRIPTOR>(iPop), 0)};
158 block.addPostProcessor(*postProcessor);
159 }
160 }
161 } // bulk indicator
162 } // iPop loop
163 setBouzidiZeroVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(block, blockGeometryStructure, iX, iY, distances);
164}

References olb::BlockLattice< T, DESCRIPTOR >::addPostProcessor(), olb::IndicatorF2D< S >::distance(), olb::BlockGeometry< T, D >::getDeltaR(), olb::BlockGeometry< T, D >::getIcGlob(), olb::BlockGeometry< T, D >::getMaterial(), olb::BlockGeometry< T, D >::getPhysR(), olb::BlockStructureD< D >::isInside(), and olb::util::sqrt().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [2/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics >
void olb::legacy::setBouzidiZeroVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 2 > & blockGeometryStructure,
int x,
int y,
int iPop,
T dist )

Definition at line 180 of file setBouzidiZeroVelocityBoundary2D.hh.

181{
182 auto postProcessor = std::unique_ptr<PostProcessorGenerator2D<T, DESCRIPTOR>>{ nullptr };
183 const Vector<int,2> c = descriptors::c<DESCRIPTOR>(iPop);
184
185 if (blockGeometryStructure.getMaterial({x-c[0], y-c[1]}) != 1) {
186 postProcessor.reset(new ZeroVelocityBounceBackPostProcessorGenerator2D<T,DESCRIPTOR>(x, y, iPop, dist));
187 }
188 else {
189 postProcessor.reset(new ZeroVelocityBouzidiLinearPostProcessorGenerator2D<T,DESCRIPTOR>(x, y, iPop, dist));
190 }
191
192 if (postProcessor && !block.isPadding({x,y})) {
193 block.addPostProcessor(*postProcessor);
194 }
195}

References olb::BlockLattice< T, DESCRIPTOR >::addPostProcessor(), olb::BlockGeometry< T, D >::getMaterial(), and olb::BlockStructureD< D >::isPadding().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [3/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics >
void olb::legacy::setBouzidiZeroVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 2 > & blockGeometryStructure,
int x,
int y,
T distances[DESCRIPTOR::q] )

Definition at line 168 of file setBouzidiZeroVelocityBoundary2D.hh.

169{
170 for (int iPop = 1; iPop < DESCRIPTOR::q ; ++iPop) {
171 if ( !util::nearZero(distances[iPop]+1) ) {
172 const Vector<int,2> c = descriptors::c<DESCRIPTOR>(iPop);
173 setBouzidiZeroVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(block, blockGeometryStructure, x-c[0], y-c[1], iPop, distances[iPop]);
174 }
175 }
176}

References olb::util::nearZero().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [4/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiZeroVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 3 > & blockGeometryStructure,
int x,
int y,
int z,
int iPop,
T dist )

Definition at line 198 of file setBouzidiZeroVelocityBoundary3D.hh.

199{
200 auto postProcessor = std::unique_ptr<PostProcessorGenerator3D<T, DESCRIPTOR>>{ nullptr };
201 const Vector<int,3> c = descriptors::c<DESCRIPTOR>(iPop);
202
203 if (blockGeometryStructure.getMaterial(x-c[0], y-c[1], z-c[2]) != 1) {
204 postProcessor.reset(new ZeroVelocityBounceBackPostProcessorGenerator3D<T,DESCRIPTOR>(x, y, z, iPop, dist));
205 }
206 else {
207 postProcessor.reset(new ZeroVelocityBouzidiLinearPostProcessorGenerator3D<T,DESCRIPTOR>(x, y, z, iPop, dist));
208 }
209
210 if (postProcessor && !block.isPadding({x,y,z})) {
211 block.addPostProcessor(*postProcessor);
212 }
213}

References olb::BlockLattice< T, DESCRIPTOR >::addPostProcessor(), olb::BlockGeometry< T, D >::getMaterial(), and olb::BlockStructureD< D >::isPadding().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [5/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiZeroVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 3 > & blockGeometryStructure,
int x,
int y,
int z,
T distances[DESCRIPTOR::q] )

Definition at line 187 of file setBouzidiZeroVelocityBoundary3D.hh.

188{
189 for (int iPop = 1; iPop < DESCRIPTOR::q; ++iPop) {
190 if ( !util::nearZero(distances[iPop]+1) ) {
191 const Vector<int,3> c = descriptors::c<DESCRIPTOR>(iPop);
192 setBouzidiZeroVelocityBoundary<T,DESCRIPTOR>(block, blockGeometryStructure, x-c[0], y-c[1], z-c[2], iPop, distances[iPop]);
193 }
194 }
195}

References olb::util::nearZero().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [6/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics >
void olb::legacy::setBouzidiZeroVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockIndicatorF2D< T > & boundaryIndicator,
BlockIndicatorF2D< T > & bulkIndicator,
IndicatorF2D< T > & geometryIndicator )

Definition at line 84 of file setBouzidiZeroVelocityBoundary2D.hh.

85{
86 block.forSpatialLocations([&](auto iX, auto iY) {
87 if (boundaryIndicator(iX,iY)) {
88 setBouzidiZeroVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(block,
89 bulkIndicator.getBlockGeometry(),
90 iX, iY,
91 bulkIndicator,
92 geometryIndicator);
93 }
94 });
95}

References olb::BlockStructureD< D >::forSpatialLocations(), and olb::BlockIndicatorF2D< T >::getBlockGeometry().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [7/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiZeroVelocityBoundary ( BlockLattice< T, DESCRIPTOR > & block,
BlockIndicatorF3D< T > & boundaryIndicator,
BlockIndicatorF3D< T > & bulkIndicator,
IndicatorF3D< T > & geometryIndicator,
T _epsFraction )

Definition at line 94 of file setBouzidiZeroVelocityBoundary3D.hh.

95{
96 block.forSpatialLocations([&](auto iX, auto iY, auto iZ) {
97 if (boundaryIndicator(iX,iY,iZ)) {
98 setBouzidiZeroVelocityBoundary1<T,DESCRIPTOR>(block,
99 boundaryIndicator.getBlockGeometry(),
100 iX, iY, iZ,
101 geometryIndicator,
102 bulkIndicator,
103 _epsFraction);
104 }
105 });
106}

References olb::BlockStructureD< D >::forSpatialLocations(), and olb::BlockIndicatorF3D< T >::getBlockGeometry().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [8/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics = BGKdynamics<T,DESCRIPTOR>>
void olb::legacy::setBouzidiZeroVelocityBoundary ( SuperLattice< T, DESCRIPTOR > & sLattice,
FunctorPtr< SuperIndicatorF2D< T > > && boundaryIndicator,
FunctorPtr< SuperIndicatorF2D< T > > && bulkIndicator,
IndicatorF2D< T > & geometryIndicator )

Initialising the BouzidiZeroVelocityBoundary on the superLattice domain.

Adds needed Cells to the Communicator _commBC in SuperLattice

Adds needed Cells to the Communicator _commBC in SuperLattice

Definition at line 53 of file setBouzidiZeroVelocityBoundary2D.hh.

56{
57 //out of superOffBoundary2D
58 T _epsFraction = 0.0001;
59 OstreamManager clout(std::cout, "setBouzidiZeroVelocityBoundary");
60 int _overlap = 1;
61
62 clout << "epsFraction=" << _epsFraction << std::endl;
63 clout.setMultiOutput(true);
64 for (int iCloc = 0; iCloc < sLattice.getLoadBalancer().size(); ++iCloc) {
65 clout << "Cuboid globiC " << sLattice.getLoadBalancer().glob(iCloc)
66 << " starts to read distances for ZeroVelocity Boundary..." << std::endl;
67 setBouzidiZeroVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice.getBlock(iCloc),
68 boundaryIndicator->getBlockIndicatorF(iCloc),
69 bulkIndicator->getBlockIndicatorF(iCloc),
70 geometryIndicator);
71 clout << "Cuboid globiC " << sLattice.getLoadBalancer().glob(iCloc)
72 << " finished reading distances for ZeroVelocity Boundary." << std::endl;
73 }
74 clout.setMultiOutput(false);
76 addPoints2CommBC<T,DESCRIPTOR>(sLattice, std::forward<decltype(boundaryIndicator)>(boundaryIndicator), _overlap);
77}

References olb::SuperLattice< T, DESCRIPTOR >::getBlock(), olb::SuperStructure< T, D >::getLoadBalancer(), and olb::OstreamManager::setMultiOutput().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [9/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiZeroVelocityBoundary ( SuperLattice< T, DESCRIPTOR > & sLattice,
FunctorPtr< SuperIndicatorF3D< T > > && boundaryIndicator,
FunctorPtr< SuperIndicatorF3D< T > > && bulkIndicator,
IndicatorF3D< T > & geometryIndicator )

Adds needed Cells to the Communicator _commBC in SuperLattice

Adds needed Cells to the Communicator _commBC in SuperLattice

Definition at line 58 of file setBouzidiZeroVelocityBoundary3D.hh.

61{
62 OstreamManager clout(std::cout, "setBouzidiZeroVelocityBoundary");
63 T _epsFraction = 0.0001;
64 int overlap = 1;
65 bool _output = false;
66 if (_output) {
67 clout << "epsFraction=" << _epsFraction << std::endl;
68 clout.setMultiOutput(true);
69 }
70 for (int iCloc = 0; iCloc < sLattice.getLoadBalancer().size(); ++iCloc) {
71 if (_output) {
72 clout << "Cuboid globiC " << sLattice.getLoadBalancer().glob(iCloc)
73 << " starts to read distances for ZeroVelocity Boundary..." << std::endl;
74 }
75 setBouzidiZeroVelocityBoundary<T,DESCRIPTOR>(sLattice.getBlock(iCloc),
76 boundaryIndicator->getBlockIndicatorF(iCloc),
77 bulkIndicator->getBlockIndicatorF(iCloc),
78 geometryIndicator, _epsFraction);
79 if (_output) {
80 clout << "Cuboid globiC " << sLattice.getLoadBalancer().glob(iCloc)
81 << " finished reading distances for ZeroVelocity Boundary." << std::endl;
82 }
83 }
84 if (_output) {
85 clout.setMultiOutput(false);
86 }
88 addPoints2CommBC(sLattice,std::forward<decltype(boundaryIndicator)>(boundaryIndicator), overlap);
89}
void addPoints2CommBC(SuperLattice< T, DESCRIPTOR > &sLattice, FunctorPtr< SuperIndicatorF2D< T > > &&indicator, int _overlap)
Adds needed Cells to the Communicator _commBC in SuperLattice.

References olb::addPoints2CommBC(), olb::SuperLattice< T, DESCRIPTOR >::getBlock(), olb::SuperStructure< T, D >::getLoadBalancer(), and olb::OstreamManager::setMultiOutput().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [10/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiZeroVelocityBoundary ( SuperLattice< T, DESCRIPTOR > & sLattice,
FunctorPtr< SuperIndicatorF3D< T > > && boundaryIndicator,
IndicatorF3D< T > & geometryIndicator,
std::vector< int > bulkMaterials )

Definition at line 48 of file setBouzidiZeroVelocityBoundary3D.hh.

50{
51 OstreamManager clout(std::cout, "setBouzidiZeroVelocityBoundary");
52 SuperGeometry<T,3>& superGeometry = boundaryIndicator->getSuperGeometry();
53 setBouzidiZeroVelocityBoundary<T,DESCRIPTOR>(sLattice, std::forward<decltype(boundaryIndicator)>(boundaryIndicator), superGeometry.getMaterialIndicator(std::move(bulkMaterials)),
54 geometryIndicator);
55}
Representation of a statistic for a parallel 2D geometry.

References olb::SuperGeometry< T, D >::getMaterialIndicator().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [11/12]

template<typename T , typename DESCRIPTOR , class MixinDynamics = BGKdynamics<T,DESCRIPTOR>>
void olb::legacy::setBouzidiZeroVelocityBoundary ( SuperLattice< T, DESCRIPTOR > & sLattice,
SuperGeometry< T, 2 > & superGeometry,
int material,
IndicatorF2D< T > & geometryIndicator,
std::vector< int > bulkMaterials = std::vector<int>(1,1) )

Set offDynamics with boundary links and post processors using indicators.

Initialising the BouzidiZeroVelocityBoundary on the superLattice domain.

Add offDynamics with initialisation of boundary links and the corresponding post processors Note: Uses information of the second neighbours of the cell (x,y) Add post processors. Ensure that offDynamics are defined!

Parameters
boundaryIndicatorIndicator describing boundary cells
bulkIndicatorIndicator describing bulk cells
geometryIndicatorIndicator describing the geometry to be bounded Initialising the BouzidiZeroVelocityBoundary on the superLattice domain

Definition at line 40 of file setBouzidiZeroVelocityBoundary2D.hh.

43{
44 setBouzidiZeroVelocityBoundary<T,DESCRIPTOR,MixinDynamics>(sLattice, superGeometry.getMaterialIndicator(material),
45 superGeometry.getMaterialIndicator(std::move(bulkMaterials)),
46 geometryIndicator);
47
48
49}

References olb::SuperGeometry< T, D >::getMaterialIndicator().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary() [12/12]

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiZeroVelocityBoundary ( SuperLattice< T, DESCRIPTOR > & sLattice,
SuperGeometry< T, 3 > & superGeometry,
int material,
IndicatorF3D< T > & geometryIndicator,
std::vector< int > bulkMaterials = std::vector<int>(1,1) )

Set offDynamics with boundary links and post processors using indicators.

Initialising the setLocalVelocityBoundary function on the superLattice domain.

Add offDynamics with initialisation of boundary links and the corresponding post processors Note: Uses information of the second neighbours of the cell (x,y,z) Add post processors. Ensure that offDynamics are defined!

Parameters
boundaryIndicatorIndicator describing boundary cells
bulkIndicatorIndicator describing bulk cells
geometryIndicatorIndicator describing the geometry to be bounded Initialising the setLocalVelocityBoundary function on the superLattice domain

Definition at line 40 of file setBouzidiZeroVelocityBoundary3D.hh.

42{
43 setBouzidiZeroVelocityBoundary<T,DESCRIPTOR>(sLattice, superGeometry.getMaterialIndicator(material),
44 geometryIndicator,
45 bulkMaterials);
46}

References olb::SuperGeometry< T, D >::getMaterialIndicator().

+ Here is the call graph for this function:

◆ setBouzidiZeroVelocityBoundary1()

template<typename T , typename DESCRIPTOR >
void olb::legacy::setBouzidiZeroVelocityBoundary1 ( BlockLattice< T, DESCRIPTOR > & block,
BlockGeometry< T, 3 > & blockGeometryStructure,
int iX,
int iY,
int iZ,
IndicatorF3D< T > & geometryIndicator,
BlockIndicatorF3D< T > & bulkIndicator,
T _epsFraction )

Definition at line 112 of file setBouzidiZeroVelocityBoundary3D.hh.

114{
115 OstreamManager clout(std::cout, "BouzidiBoundary");
116 T distances[DESCRIPTOR::q];
117 std::fill(std::begin(distances), std::end(distances), -1);
118
119 for (int iPop = 1; iPop < DESCRIPTOR::q ; ++iPop) {
120 const Vector<int,3> c = descriptors::c<DESCRIPTOR>(iPop);
121 const int iXn = iX + c[0];
122 const int iYn = iY + c[1];
123 const int iZn = iZ + c[2];
124 if (blockGeometryStructure.isInside(iXn,iYn,iZn)) {
125 if (bulkIndicator(iXn,iYn,iZn)) {
126 T dist = -1;
127 T physR[3];
128 blockGeometryStructure.getPhysR(physR,{iXn,iYn,iZn});
129 T voxelSize=blockGeometryStructure.getDeltaR();
130
131 Vector<T,3> physC(physR);
132
133 Vector<T,3> direction(-voxelSize*c[0],-voxelSize*c[1],-voxelSize*c[2]);
134 T cPhysNorm = voxelSize*util::sqrt(c[0]*c[0]+c[1]*c[1]+c[2]*c[2]);
135
136 if (!geometryIndicator.distance(dist,physC,direction,blockGeometryStructure.getIcGlob() ) ) {
137 T epsX = voxelSize*c[0]*_epsFraction;
138 T epsY = voxelSize*c[1]*_epsFraction;
139 T epsZ = voxelSize*c[2]*_epsFraction;
140
141 Vector<T,3> physC2(physC);
142 physC2[0] += epsX;
143 physC2[1] += epsY;
144 physC2[2] += epsZ;
145 Vector<T,3> direction2(direction);
146 direction2[0] -= 2.*epsX;
147 direction2[1] -= 2.*epsY;
148 direction2[2] -= 2.*epsZ;
149
150 if ( !geometryIndicator.distance(dist,physC2,direction2,blockGeometryStructure.getIcGlob())) {
151 clout << "ERROR: no boundary found at (" << iXn << "," << iYn << "," << iZn <<") ~ ("
152 << physR[0] << "," << physR[1] << "," << physR[2] <<"), "
153 << "in direction " << descriptors::opposite<DESCRIPTOR >(iPop)
154 << std::endl;
155 break;
156 }
157 T distNew = (dist - util::sqrt(epsX*epsX+epsY*epsY+epsZ*epsZ))/cPhysNorm;
158 if (distNew < 0.5) {
159 dist = 0;
160 }
161 else {
162 dist = 0.5 * cPhysNorm;
163 clout << "WARNING: distance at (" << iXn << "," << iYn << "," << iZn <<") ~ ("
164 << physR[0] << "," << physR[1] << "," << physR[2] <<"), "
165 << "in direction " << descriptors::opposite<DESCRIPTOR >(iPop) << ": "
166 << distNew
167 << " rounded to "
168 << dist/cPhysNorm
169 << std::endl;
170 }
171 }
172 distances[descriptors::opposite<DESCRIPTOR>(iPop)] = dist/cPhysNorm;
173 }
174 else {
175 if (blockGeometryStructure.getMaterial(iXn,iYn,iZn) != 0) {
176 auto postProcessor = std::unique_ptr<PostProcessorGenerator3D<T, DESCRIPTOR>>{
177 new ZeroVelocityBounceBackPostProcessorGenerator3D<T,DESCRIPTOR>(iXn, iYn, iZn, descriptors::opposite<DESCRIPTOR>(iPop), 0) };
178 block.addPostProcessor(*postProcessor);
179 }
180 }
181 }
182 }
183 setBouzidiZeroVelocityBoundary<T,DESCRIPTOR>(block,blockGeometryStructure, iX, iY, iZ, distances);
184}

References olb::BlockLattice< T, DESCRIPTOR >::addPostProcessor(), olb::IndicatorF3D< S >::distance(), olb::BlockGeometry< T, D >::getDeltaR(), olb::BlockGeometry< T, D >::getIcGlob(), olb::BlockGeometry< T, D >::getMaterial(), olb::BlockGeometry< T, D >::getPhysR(), olb::BlockStructureD< D >::isInside(), and olb::util::sqrt().

+ Here is the call graph for this function:

◆ setOffDynamics()

template<typename T , typename DESCRIPTOR , class MixinDynamics >
void olb::legacy::setOffDynamics ( BlockLattice< T, DESCRIPTOR > & block,
int x,
int y,
T location[DESCRIPTOR::d],
T distances[DESCRIPTOR::q] )

Definition at line 210 of file setBouzidiVelocityBoundary2D.hh.

211{
212 //Dynamics<T,DESCRIPTOR>* dynamics = BoundaryManager::getOffDynamics(location, distances);
213 Dynamics<T,DESCRIPTOR>* dynamics = new legacy::OffDynamics<T, DESCRIPTOR>(location, distances);
214 block.defineDynamics({x,y}, dynamics);
215}
Interface for per-cell dynamics.
Definition interface.h:54

References olb::BlockLattice< T, DESCRIPTOR >::defineDynamics().

+ Here is the call graph for this function: