OpenLB 1.7
Loading...
Searching...
No Matches
Public Types | Public Member Functions | List of all members
olb::collision::KrauseH< COLLISION >::type< DESCRIPTOR, MOMENTA, EQUILIBRIUM > Struct Template Reference

#include <porousBGKdynamics.h>

+ Collaboration diagram for olb::collision::KrauseH< COLLISION >::type< DESCRIPTOR, MOMENTA, EQUILIBRIUM >:

Public Types

using MomentaF = typename MOMENTA::template type<DESCRIPTOR>
 
using EquilibriumF = typename EQUILIBRIUM::template type<DESCRIPTOR,MOMENTA>
 
using CollisionO = typename COLLISION::template type<DESCRIPTOR,MOMENTA,EQUILIBRIUM>
 

Public Member Functions

template<typename CELL , typename PARAMETERS , typename V = typename CELL::value_t>
CellStatistic< V > apply (CELL &cell, PARAMETERS &parameters)
 

Detailed Description

template<typename COLLISION>
template<typename DESCRIPTOR, typename MOMENTA, typename EQUILIBRIUM>
struct olb::collision::KrauseH< COLLISION >::type< DESCRIPTOR, MOMENTA, EQUILIBRIUM >

Definition at line 346 of file porousBGKdynamics.h.

Member Typedef Documentation

◆ CollisionO

template<typename COLLISION >
template<typename DESCRIPTOR , typename MOMENTA , typename EQUILIBRIUM >
using olb::collision::KrauseH< COLLISION >::type< DESCRIPTOR, MOMENTA, EQUILIBRIUM >::CollisionO = typename COLLISION::template type<DESCRIPTOR,MOMENTA,EQUILIBRIUM>

Definition at line 349 of file porousBGKdynamics.h.

◆ EquilibriumF

template<typename COLLISION >
template<typename DESCRIPTOR , typename MOMENTA , typename EQUILIBRIUM >
using olb::collision::KrauseH< COLLISION >::type< DESCRIPTOR, MOMENTA, EQUILIBRIUM >::EquilibriumF = typename EQUILIBRIUM::template type<DESCRIPTOR,MOMENTA>

Definition at line 348 of file porousBGKdynamics.h.

◆ MomentaF

template<typename COLLISION >
template<typename DESCRIPTOR , typename MOMENTA , typename EQUILIBRIUM >
using olb::collision::KrauseH< COLLISION >::type< DESCRIPTOR, MOMENTA, EQUILIBRIUM >::MomentaF = typename MOMENTA::template type<DESCRIPTOR>

Definition at line 347 of file porousBGKdynamics.h.

Member Function Documentation

◆ apply()

template<typename COLLISION >
template<typename DESCRIPTOR , typename MOMENTA , typename EQUILIBRIUM >
template<typename CELL , typename PARAMETERS , typename V = typename CELL::value_t>
CellStatistic< V > olb::collision::KrauseH< COLLISION >::type< DESCRIPTOR, MOMENTA, EQUILIBRIUM >::apply ( CELL & cell,
PARAMETERS & parameters )
inline

Molecular realaxation time

Turbulent realaxation time

Effective realaxation time

Definition at line 352 of file porousBGKdynamics.h.

352 {
353 V rho, u[DESCRIPTOR::d];
354 V newOmega[DESCRIPTOR::d];
355 V _fieldTmp[4];
356 _fieldTmp[0] = V{1};
357 _fieldTmp[1] = V{};
358 _fieldTmp[2] = V{};
359 _fieldTmp[3] = V{};
360
361 const V smagoConst = parameters.template get<collision::LES::Smagorinsky>();
362 V preFactor = smagoConst*smagoConst
363 * descriptors::invCs2<V,DESCRIPTOR>()*descriptors::invCs2<V,DESCRIPTOR>()
364 * 2*util::sqrt(2);
365
366 MomentaF().computeRhoU(cell, rho, u);
367
368 // compute newOmega
369 V omega = parameters.template get<descriptors::OMEGA>();
370 V uSqr = u[0]*u[0];
371 for (int iDim=0; iDim<DESCRIPTOR::d; iDim++) {
372 uSqr += u[iDim]*u[iDim];
373 }
375 V tau_mol = 1./omega;
376 for (int iPop=0; iPop<DESCRIPTOR::q; iPop++) {
377 V fNeq = util::fabs(cell[iPop] - EquilibriumF().compute(iPop, rho, u));
379 V tau_turb = 0.5*(util::sqrt(tau_mol*tau_mol+(preFactor*fNeq))-tau_mol);
381 V tau_eff = tau_mol + tau_turb;
382 newOmega[iPop] = 1./tau_eff;
383 }
384 parameters.template set<descriptors::OMEGA>(newOmega);
385
386
387 V vel_denom = cell.template getField<descriptors::VELOCITY_DENOMINATOR>();
388 if (vel_denom > std::numeric_limits<V>::epsilon()) {
389 V porosity = cell.template getField<descriptors::POROSITY>(); // prod(1-smoothInd)
390 auto vel_num = cell.template getFieldPointer<descriptors::VELOCITY_NUMERATOR>();
391 porosity = 1.-porosity; // 1-prod(1-smoothInd)
392 for (int i=0; i<DESCRIPTOR::d; i++) {
393 u[i] += porosity * (vel_num[i] / vel_denom - u[i]);
394 }
395 }
396 uSqr = CollisionO().apply(cell, parameters);
397
398 cell.template setField<descriptors::POROSITY>(_fieldTmp[0]);
399 cell.template setField<descriptors::VELOCITY_NUMERATOR>({_fieldTmp[1], _fieldTmp[2]});
400 cell.template setField<descriptors::VELOCITY_DENOMINATOR>(_fieldTmp[3]);
401
402 return {rho, uSqr};
403 }
cpu::simd::Pack< T > sqrt(cpu::simd::Pack< T > value)
Definition pack.h:100
cpu::simd::Pack< T > fabs(cpu::simd::Pack< T > value)
Definition pack.h:106
typename COLLISION::template type< DESCRIPTOR, MOMENTA, EQUILIBRIUM > CollisionO
typename EQUILIBRIUM::template type< DESCRIPTOR, MOMENTA > EquilibriumF
typename MOMENTA::template type< DESCRIPTOR > MomentaF
meta::list< descriptors::OMEGA, LES::Smagorinsky > parameters

References olb::util::fabs(), and olb::util::sqrt().

+ Here is the call graph for this function:

The documentation for this struct was generated from the following file: