oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM > Class Template Reference

#include <generalised_newtonian_constitutive_models.h>

+ Inheritance diagram for oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >:

Public Member Functions

 HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation (double *yield_stress_pt, double *flow_index_pt, double *critical_second_invariant_pt)
 "Cutoff regularised" Herschel Bulkley constitutive equation More...
 
double calculate_cutoff_viscosity ()
 
double calculate_zero_shear_viscosity ()
 Function that calculates the viscosity at zero I2. More...
 
void report_cut_off_values (double &cut_off_invariant, double &cut_off_viscosity, double &zero_shear_viscosity)
 Report cutoff values. More...
 
double viscosity (const double &second_invariant_of_rate_of_strain_tensor)
 Viscosity ratio as a fct of strain rate invariant. More...
 
double dviscosity_dinvariant (const double &second_invariant_of_rate_of_strain_tensor)
 Deriv of viscosity w.r.t. strain rate invariant. More...
 
- Public Member Functions inherited from oomph::GeneralisedNewtonianConstitutiveEquation< DIM >
 GeneralisedNewtonianConstitutiveEquation ()
 Empty constructor. More...
 
virtual ~GeneralisedNewtonianConstitutiveEquation ()
 Empty virtual destructor. More...
 

Private Attributes

doubleYield_stress_pt
 yield stress tau_y More...
 
doubleFlow_index_pt
 power law index n More...
 
doubleCritical_second_invariant_pt
 
double a
 
double b
 
double c
 
double alpha
 

Detailed Description

template<unsigned DIM>
class oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >

A GeneralisedNewtonianConstitutiveEquation class defining a Herschel-Bulkley fluid using Tanner and Milthorpe's (1983) regularisation with a smooth transition using a quadratic

Constructor & Destructor Documentation

◆ HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation()

template<unsigned DIM>
oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation ( double yield_stress_pt,
double flow_index_pt,
double critical_second_invariant_pt 
)
inline

"Cutoff regularised" Herschel Bulkley constitutive equation

Blend over one order of magnitude

Calculate the coefficients of the quadratic equation

get the cutoff viscosity

get the zero shear viscosity

392  : GeneralisedNewtonianConstitutiveEquation<DIM>(),
393  Yield_stress_pt(yield_stress_pt),
394  Flow_index_pt(flow_index_pt),
395  Critical_second_invariant_pt(critical_second_invariant_pt),
396  a(0.0),
397  b(0.0),
398  c(0.0)
399  {
401  alpha = 0.1;
402 
405  {
406  a = (pow(2.0, *Flow_index_pt - 3.0) * (*Flow_index_pt - 1.0) *
408  (*Flow_index_pt - 1.0) / 2.0 - 2.0) -
409  (*Yield_stress_pt) /
410  (8.0 * pow(fabs(*Critical_second_invariant_pt), 5.0 / 2.0))) /
411  (1.0 - alpha);
412  b = -2.0 * a * alpha * fabs(*Critical_second_invariant_pt);
413  c = (*Yield_stress_pt) /
416  *Flow_index_pt - 1.0) -
419  }
420 
422  double cut_off_viscosity = calculate_cutoff_viscosity();
423 
425  double zero_shear_viscosity = calculate_zero_shear_viscosity();
426 
427  oomph_info << "HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation: "
428  << " zero shear viscosity = " << zero_shear_viscosity
429  << std::endl;
430 
431  oomph_info << "HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation: "
432  << " cut off viscosity = " << cut_off_viscosity << std::endl;
433 
434  oomph_info << " "
435  << " cutoff invariant = " << *Critical_second_invariant_pt
436  << std::endl;
437  }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
double * Flow_index_pt
power law index n
Definition: generalised_newtonian_constitutive_models.h:366
double * Critical_second_invariant_pt
Definition: generalised_newtonian_constitutive_models.h:370
double calculate_cutoff_viscosity()
Definition: generalised_newtonian_constitutive_models.h:441
double a
Definition: generalised_newtonian_constitutive_models.h:377
double alpha
Definition: generalised_newtonian_constitutive_models.h:383
double calculate_zero_shear_viscosity()
Function that calculates the viscosity at zero I2.
Definition: generalised_newtonian_constitutive_models.h:454
double b
Definition: generalised_newtonian_constitutive_models.h:378
double c
Definition: generalised_newtonian_constitutive_models.h:379
double * Yield_stress_pt
yield stress tau_y
Definition: generalised_newtonian_constitutive_models.h:363
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::a, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::alpha, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::b, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::c, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::calculate_cutoff_viscosity(), oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::calculate_zero_shear_viscosity(), oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::Critical_second_invariant_pt, boost::multiprecision::fabs(), oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::Flow_index_pt, oomph::oomph_info, Eigen::bfloat16_impl::pow(), and sqrt().

Member Function Documentation

◆ calculate_cutoff_viscosity()

◆ calculate_zero_shear_viscosity()

◆ dviscosity_dinvariant()

template<unsigned DIM>
double oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::dviscosity_dinvariant ( const double second_invariant_of_rate_of_strain_tensor)
inlinevirtual

Deriv of viscosity w.r.t. strain rate invariant.

Implements oomph::GeneralisedNewtonianConstitutiveEquation< DIM >.

498  {
499  if (fabs(second_invariant_of_rate_of_strain_tensor) <
501  {
502  return 0.0;
503  }
504  else if (fabs(second_invariant_of_rate_of_strain_tensor) <
506  {
507  return 2.0 * a * fabs(second_invariant_of_rate_of_strain_tensor) + b;
508  }
509 
510  return pow(2.0, *Flow_index_pt - 2.0) * (*Flow_index_pt - 1.0) *
511  pow(fabs(second_invariant_of_rate_of_strain_tensor),
512  (*Flow_index_pt - 1.0) / 2.0 - 1.0) -
513  (*Yield_stress_pt) /
514  (4.0 * pow(fabs(second_invariant_of_rate_of_strain_tensor),
515  3.0 / 2.0));
516  }

References oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::a, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::alpha, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::b, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::Critical_second_invariant_pt, boost::multiprecision::fabs(), oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::Flow_index_pt, and Eigen::bfloat16_impl::pow().

◆ report_cut_off_values()

◆ viscosity()

template<unsigned DIM>
double oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::viscosity ( const double second_invariant_of_rate_of_strain_tensor)
inlinevirtual

Viscosity ratio as a fct of strain rate invariant.

Implements oomph::GeneralisedNewtonianConstitutiveEquation< DIM >.

476  {
477  if (fabs(second_invariant_of_rate_of_strain_tensor) <
479  {
481  }
482  else if (fabs(second_invariant_of_rate_of_strain_tensor) <
484  {
485  return a * pow(fabs(second_invariant_of_rate_of_strain_tensor), 2.0) +
486  b * fabs(second_invariant_of_rate_of_strain_tensor) + c;
487  }
488 
489  return (*Yield_stress_pt) /
490  (2.0 * sqrt(fabs(second_invariant_of_rate_of_strain_tensor))) +
491  pow(2.0 * sqrt(fabs(second_invariant_of_rate_of_strain_tensor)),
492  *Flow_index_pt - 1.0);
493  }

References oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::a, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::alpha, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::b, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::c, oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::calculate_zero_shear_viscosity(), oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::Critical_second_invariant_pt, boost::multiprecision::fabs(), oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::Flow_index_pt, Eigen::bfloat16_impl::pow(), sqrt(), and oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::Yield_stress_pt.

Member Data Documentation

◆ a

We use a quadratic function to smoothly blend from the Herschel Bulkley model at the cut-off to a constant viscosity as the strain rate /approaches zero; this way we avoid the discontinuity of the gradient at the cut-off, which is present in the classic Tanner Milthorpe regularisation

Referenced by oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::calculate_zero_shear_viscosity(), oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::dviscosity_dinvariant(), oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation(), and oomph::HerschelBulkleyTanMilRegWithBlendingConstitutiveEquation< DIM >::viscosity().

◆ alpha

◆ b

◆ c

◆ Critical_second_invariant_pt

◆ Flow_index_pt

◆ Yield_stress_pt


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