OpenLB 1.7
Loading...
Searching...
No Matches
Public Member Functions | List of all members
olb::SuperLatticeYplus3D< T, DESCRIPTOR > Class Template Reference

functor to get pointwise yPlus from rho, shear stress and local density on local lattices More...

#include <turbulentF3D.h>

+ Inheritance diagram for olb::SuperLatticeYplus3D< T, DESCRIPTOR >:
+ Collaboration diagram for olb::SuperLatticeYplus3D< T, DESCRIPTOR >:

Public Member Functions

 SuperLatticeYplus3D (SuperLattice< T, DESCRIPTOR > &sLattice, const UnitConverter< T, DESCRIPTOR > &converter, SuperGeometry< T, 3 > &superGeometry, IndicatorF3D< T > &indicator, const int material)
 
bool operator() (T output[], const int input[]) override
 
- Public Member Functions inherited from olb::SuperLatticePhysF3D< T, DESCRIPTOR >
UnitConverter< T, DESCRIPTOR > const & getConverter () const
 
- Public Member Functions inherited from olb::SuperLatticeF3D< T, DESCRIPTOR >
SuperLattice< T, DESCRIPTOR > & getSuperLattice ()
 
bool operator() (T output[], const int input[])
 
- Public Member Functions inherited from olb::SuperF3D< T, T >
SuperF3D< T, T > & operator- (SuperF3D< T, T > &rhs)
 
SuperF3D< T, T > & operator+ (SuperF3D< T, T > &rhs)
 
SuperF3D< T, T > & operator* (SuperF3D< T, T > &rhs)
 
SuperF3D< T, T > & operator/ (SuperF3D< T, T > &rhs)
 
SuperStructure< T, 3 > & getSuperStructure ()
 
int getBlockFSize () const
 
BlockF3D< T > & getBlockF (int iCloc)
 
bool operator() (T output[], const int input[])
 
- Public Member Functions inherited from olb::GenericF< T, S >
virtual ~GenericF ()=default
 
int getSourceDim () const
 read only access to member variable _m
 
int getTargetDim () const
 read only access to member variable _n
 
std::string & getName ()
 read and write access to name
 
std::string const & getName () const
 read only access to name
 
virtual bool operator() (T output[], const S input[])=0
 has to be implemented for 'every' derived class
 
bool operator() (T output[])
 wrapper that call the pure virtual operator() (T output[], const S input[]) from above
 
bool operator() (T output[], S input0)
 
bool operator() (T output[], S input0, S input1)
 
bool operator() (T output[], S input0, S input1, S input2)
 
bool operator() (T output[], S input0, S input1, S input2, S input3)
 

Additional Inherited Members

- Public Types inherited from olb::SuperLatticeF3D< T, DESCRIPTOR >
using identity_functor_type = SuperLatticeIdentity3D<T,DESCRIPTOR>
 
- Public Types inherited from olb::SuperF3D< T, T >
using identity_functor_type
 
- Public Types inherited from olb::GenericF< T, S >
using targetType = T
 
using sourceType = S
 
- Public Attributes inherited from olb::GenericF< T, S >
std::shared_ptr< GenericF< T, S > > _ptrCalcC
 memory management, frees resouces (calcClass)
 
- Static Public Attributes inherited from olb::SuperF3D< T, T >
static constexpr bool isSuper
 
static constexpr unsigned d
 
- Protected Member Functions inherited from olb::SuperLatticePhysF3D< T, DESCRIPTOR >
 SuperLatticePhysF3D (SuperLattice< T, DESCRIPTOR > &sLattice, const UnitConverter< T, DESCRIPTOR > &converter, int targetDim)
 
- Protected Member Functions inherited from olb::SuperLatticeF3D< T, DESCRIPTOR >
 SuperLatticeF3D (SuperLattice< T, DESCRIPTOR > &superLattice, int targetDim)
 
- Protected Member Functions inherited from olb::SuperF3D< T, T >
 SuperF3D (SuperStructure< T, 3 > &superStructure, int targetDim)
 
- Protected Member Functions inherited from olb::GenericF< T, S >
 GenericF (int targetDim, int sourceDim)
 
- Protected Attributes inherited from olb::SuperLatticePhysF3D< T, DESCRIPTOR >
const UnitConverter< T, DESCRIPTOR > & _converter
 
- Protected Attributes inherited from olb::SuperLatticeF3D< T, DESCRIPTOR >
SuperLattice< T, DESCRIPTOR > & _sLattice
 
- Protected Attributes inherited from olb::SuperF3D< T, T >
SuperStructure< T, 3 > & _superStructure
 
std::vector< std::unique_ptr< BlockF3D< T > > > _blockF
 Super functors may consist of several BlockF3D<W> derived functors.
 

Detailed Description

template<typename T, typename DESCRIPTOR>
class olb::SuperLatticeYplus3D< T, DESCRIPTOR >

functor to get pointwise yPlus from rho, shear stress and local density on local lattices

Definition at line 47 of file turbulentF3D.h.

Constructor & Destructor Documentation

◆ SuperLatticeYplus3D()

template<typename T , typename DESCRIPTOR >
olb::SuperLatticeYplus3D< T, DESCRIPTOR >::SuperLatticeYplus3D ( SuperLattice< T, DESCRIPTOR > & sLattice,
const UnitConverter< T, DESCRIPTOR > & converter,
SuperGeometry< T, 3 > & superGeometry,
IndicatorF3D< T > & indicator,
const int material )

Definition at line 34 of file turbulentF3D.hh.

37 : SuperLatticePhysF3D<T,DESCRIPTOR>(sLattice,converter,1),
38 _superGeometry(superGeometry), _indicator(indicator), _material(material)
39{
40 this->getName() = "yPlus";
41}
std::string & getName()
read and write access to name
Definition genericF.hh:51

References olb::GenericF< T, S >::getName().

+ Here is the call graph for this function:

Member Function Documentation

◆ operator()()

template<typename T , typename DESCRIPTOR >
bool olb::SuperLatticeYplus3D< T, DESCRIPTOR >::operator() ( T output[],
const int input[] )
override

Definition at line 44 of file turbulentF3D.hh.

45{
46 int globIC = input[0];
47 int locix = input[1];
48 int lociy = input[2];
49 int lociz = input[3];
50
51 output[0]=T();
52 if ( this->_sLattice.getLoadBalancer().rank(globIC) == singleton::mpi().getRank() ) {
53 std::vector<T> normalTemp(3,T());
54 std::vector<T> normal(3,T()); // normalized
55 unsigned counter {0};
56 T distance = T();
57 if (_superGeometry.get(input) == 1) {
58 for (int iPop = 1; iPop < DESCRIPTOR::q; ++iPop) {
59 if (_superGeometry.get(input[0],
60 input[1] + descriptors::c<DESCRIPTOR>(iPop,0),
61 input[2] + descriptors::c<DESCRIPTOR>(iPop,1),
62 input[3] + descriptors::c<DESCRIPTOR>(iPop,2)) == _material) {
63 ++counter;
64 normalTemp[0] += descriptors::c<DESCRIPTOR>(iPop,0);
65 normalTemp[1] += descriptors::c<DESCRIPTOR>(iPop,1);
66 normalTemp[2] += descriptors::c<DESCRIPTOR>(iPop,2);
67 }
68 }
69 if ( counter > 0 ) {
70 // get physical Coordinates at intersection
71
72 std::vector<T> physR(3, T());
73 _superGeometry.getCuboidGeometry().getPhysR(&(physR[0]), &(input[0]));
74
75 T voxelSize = _superGeometry.getCuboidGeometry().get(globIC).getDeltaR();
76
77 normal = util::normalize(normalTemp);
78
79 std::vector<T> direction(normal);
80 direction[0] = voxelSize*normal[0]*2.;
81 direction[1] = voxelSize*normal[1]*2.;
82 direction[2] = voxelSize*normal[2]*2.;
83
84 // calculate distance to STL file
85 if ( _indicator.distance(distance, physR, direction) ) {
86 // call stress at this point
87 T rho;
88 T u[3];
89 T pi[6];
90 auto cell = this->_sLattice.getBlock(this->_sLattice.getLoadBalancer().loc(globIC)).get(locix, lociy, lociz);
91 cell.computeAllMomenta(rho, u, pi);
92
93 // Compute phys stress tau = mu*du/dx
94 T omega = 1. / this->_converter.getLatticeRelaxationTime();
95 T dt = this->_converter.getConversionFactorTime();
96 T physFactor = -omega*descriptors::invCs2<T,DESCRIPTOR>()/rho/2./dt*this->_converter.getPhysDensity(rho)*this->_converter.getPhysViscosity();
97
98 // Totel Stress projected from cell in normal direction on obstacle
99 T Rx = pi[0]*physFactor*normal[0] + pi[1]*physFactor*normal[1] + pi[2]*physFactor*normal[2];
100 T Ry = pi[1]*physFactor*normal[0] + pi[3]*physFactor*normal[1] + pi[4]*physFactor*normal[2];
101 T Rz = pi[2]*physFactor*normal[0] + pi[4]*physFactor*normal[1] + pi[5]*physFactor*normal[2];
102
103 // Stress appearing as pressure in corresponding direction is calculated and substracted
104 T R_res_pressure = normal[0]*pi[0]*physFactor*normal[0] + normal[0]*pi[1]*physFactor*normal[1] + normal[0]*pi[2]*physFactor*normal[2]
105 +normal[1]*pi[1]*physFactor*normal[0] + normal[1]*pi[3]*physFactor*normal[1] + normal[1]*pi[4]*physFactor*normal[2]
106 +normal[2]*pi[2]*physFactor*normal[0] + normal[2]*pi[4]*physFactor*normal[1] + normal[2]*pi[5]*physFactor*normal[2];
107
108 Rx -= R_res_pressure *normal[0];
109 Ry -= R_res_pressure *normal[1];
110 Rz -= R_res_pressure *normal[2];
111
112 T tau_wall = util::sqrt(Rx*Rx+Ry*Ry+Rz*Rz);
113 T u_tau = util::sqrt(tau_wall/this->_converter.getPhysDensity(rho));
114 //y_plus
115 output[0] = distance*u_tau / this->_converter.getPhysViscosity();
116 } // if 4
117 }
118 }
119 }
120 return true;
121}
virtual bool distance(S &distance, const Vector< S, 3 > &origin, S precision, const Vector< S, 3 > &direction)
int get(int iCglob, LatticeR< D > latticeR) const
Read only access to the material numbers, error handling: returns 0 if data is not available.
SuperLattice< T, DESCRIPTOR > & _sLattice
const UnitConverter< T, DESCRIPTOR > & _converter
CuboidGeometry< T, D > & getCuboidGeometry()
Read and write access to cuboid geometry.
int getRank() const
Returns the process ID.
MpiManager & mpi()
cpu::simd::Pack< T > sqrt(cpu::simd::Pack< T > value)
Definition pack.h:100
bool distance(S &distance, const Vector< S, D > &origin, const Vector< S, D > &direction, S precision, S pitch, F1 isInside, F2 isInsideBoundingBox)
Vector< T, D > normalize(const Vector< T, D > &a)

References olb::singleton::MpiManager::getRank(), olb::singleton::mpi(), olb::util::normalize(), and olb::util::sqrt().

+ Here is the call graph for this function:

The documentation for this class was generated from the following files: