HeatFluidCoupledSpecies< NormalForceSpecies > Class Template Reference

Species for the HeatFluidCoupledParticle. More...

#include <HeatFluidCoupledSpecies.h>

+ Inheritance diagram for HeatFluidCoupledSpecies< NormalForceSpecies >:

Public Types

typedef HeatFluidCoupledInteraction< typename NormalForceSpecies::InteractionType > InteractionType
 
- Public Types inherited from ThermalSpecies< NormalForceSpecies >
typedef ThermalInteraction< typename NormalForceSpecies::InteractionType > InteractionType
 The correct Interaction type for this FrictionForceSpecies. More...
 

Public Member Functions

 HeatFluidCoupledSpecies ()
 The default constructor. More...
 
 HeatFluidCoupledSpecies (const HeatFluidCoupledSpecies &s)
 The default copy constructor. More...
 
virtual ~HeatFluidCoupledSpecies ()
 The default destructor. More...
 
void write (std::ostream &os) const
 Writes the species properties to an output stream. More...
 
void read (std::istream &is)
 Reads the species properties from an input stream. More...
 
std::string getBaseName () const
 Used in Species::getName to obtain a unique name for each Species. More...
 
Mdouble getMassTransferCoefficient () const
 Allows massTransferCoefficient_ to be accessed. More...
 
void setMassTransferCoefficient (Mdouble massTransferCoefficient)
 Allows massTransferCoefficient_ to be changed. More...
 
Mdouble getLatentHeatVaporization () const
 Allows latentHeatVaporization_ to be accessed. More...
 
void setLatentHeatVaporization (Mdouble latentHeatVaporization)
 Allows latentHeatVaporization_ to be changed. More...
 
Mdouble getLiquidDensity () const
 Allows liquidDensity_ to be accessed. More...
 
void setLiquidDensity (Mdouble liquidDensity)
 Allows liquidDensity_ to be changed. More...
 
Mdouble getEvaporationCoefficientA () const
 Allows evaporationCoefficientA_ to be accessed. More...
 
void setEvaporationCoefficientA (Mdouble evaporationCoefficientA)
 Allows evaporationCoefficientA_ to be changed. More...
 
Mdouble getEvaporationCoefficientB () const
 Allows evaporationCoefficientB_ to be accessed. More...
 
void setEvaporationCoefficientB (Mdouble evaporationCoefficientB)
 Allows evaporationCoefficientB_ to be changed. More...
 
Mdouble getAmbientHumidity () const
 Allows ambientHumidity_ to be accessed. More...
 
void setAmbientHumidity (Mdouble ambientHumidity)
 Allows ambientHumidity_ to be changed. More...
 
Mdouble getAmbientEquilibriumMoistureContent () const
 Allows ambientEquilibriumMoistureContent_ to be accessed. More...
 
void setAmbientEquilibriumMoistureContent (Mdouble ambientEquilibriumMoistureContent)
 Allows ambientEquilibriumMoistureContent_ to be changed. More...
 
Mdouble getAmbientVapourConcentration () const
 Allows ambientVapourConcentration_ to be accessed. More...
 
void setAmbientVapourConcentration (Mdouble ambientVapourConcentration)
 Allows ambientVapourConcentration_ to be changed. More...
 
Mdouble getAmbientTemperature () const
 Allows ambientTemperature_ to be accessed. More...
 
void setAmbientTemperature (Mdouble ambientTemperature)
 Allows ambientTemperature_ to be changed. More...
 
void actionsAfterTimeStep (BaseParticle *particle) const override
 
std::array< double, 2 > f (double liquidVolume, double temperature, double mass, double surfaceArea) const
 f1 is used in Runge–Kutta method. More...
 
- Public Member Functions inherited from ThermalSpecies< NormalForceSpecies >
 ThermalSpecies ()
 The default constructor. More...
 
 ThermalSpecies (const ThermalSpecies &s)
 The default copy constructor. More...
 
virtual ~ThermalSpecies ()
 The default destructor. More...
 
void write (std::ostream &os) const
 Writes the species properties to an output stream. More...
 
void read (std::istream &is)
 Reads the species properties from an input stream. More...
 
std::string getBaseName () const
 Used in Species::getName to obtain a unique name for each Species. More...
 
Mdouble getHeatCapacity () const
 Allows heatCapacity_ to be accessed. More...
 
void setHeatCapacity (Mdouble heatCapacity)
 Allows heatCapacity_ to be changed. More...
 
Mdouble getThermalConductivity () const
 Allows heatCapacity_ to be accessed. More...
 
void setThermalConductivity (Mdouble thermalConductivity)
 Allows heatCapacity_ to be changed. More...
 

Private Attributes

Mdouble massTransferCoefficient_
 The mass transfer rate (m/s) More...
 
Mdouble latentHeatVaporization_
 The latent heat of vaporization (J/kg) More...
 
Mdouble liquidDensity_
 The liquid density (kg/m^3) More...
 
Mdouble evaporationCoefficientA_
 The evaporation coefficient a (dimensionless) More...
 
Mdouble evaporationCoefficientB_
 The evaporation coefficient b (dimensionless) More...
 
Mdouble ambientHumidity_
 The ambient humidity (dimensionless, between 0 and 1, but cannot be 0) More...
 
Mdouble ambientEquilibriumMoistureContent_
 The ambient equilibrium moisture content (dimensionless, between 0 and 1). More...
 
Mdouble ambientVapourConcentration_
 The ambient vapour concentration (kg/m^3). More...
 
Mdouble ambientTemperature_
 The ambient temperature (K). More...
 

Detailed Description

template<class NormalForceSpecies>
class HeatFluidCoupledSpecies< NormalForceSpecies >

Species for the HeatFluidCoupledParticle.

Member Typedef Documentation

◆ InteractionType

template<class NormalForceSpecies >
typedef HeatFluidCoupledInteraction<typename NormalForceSpecies::InteractionType> HeatFluidCoupledSpecies< NormalForceSpecies >::InteractionType

Constructor & Destructor Documentation

◆ HeatFluidCoupledSpecies() [1/2]

template<class NormalForceSpecies >
HeatFluidCoupledSpecies< NormalForceSpecies >::HeatFluidCoupledSpecies

The default constructor.

169 {
172  liquidDensity_=1.0;
175  // note, this default value is unrealistically high; 0.3 is realistic.
176  ambientHumidity_=1.0;
180 }
Mdouble ambientVapourConcentration_
The ambient vapour concentration (kg/m^3).
Definition: HeatFluidCoupledSpecies.h:157
Mdouble evaporationCoefficientB_
The evaporation coefficient b (dimensionless)
Definition: HeatFluidCoupledSpecies.h:142
Mdouble massTransferCoefficient_
The mass transfer rate (m/s)
Definition: HeatFluidCoupledSpecies.h:122
Mdouble ambientHumidity_
The ambient humidity (dimensionless, between 0 and 1, but cannot be 0)
Definition: HeatFluidCoupledSpecies.h:147
Mdouble ambientEquilibriumMoistureContent_
The ambient equilibrium moisture content (dimensionless, between 0 and 1).
Definition: HeatFluidCoupledSpecies.h:152
Mdouble evaporationCoefficientA_
The evaporation coefficient a (dimensionless)
Definition: HeatFluidCoupledSpecies.h:137
Mdouble latentHeatVaporization_
The latent heat of vaporization (J/kg)
Definition: HeatFluidCoupledSpecies.h:127
Mdouble liquidDensity_
The liquid density (kg/m^3)
Definition: HeatFluidCoupledSpecies.h:132
Mdouble ambientTemperature_
The ambient temperature (K).
Definition: HeatFluidCoupledSpecies.h:162
Definition: ThermalSpecies.h:35

References HeatFluidCoupledSpecies< NormalForceSpecies >::ambientEquilibriumMoistureContent_, HeatFluidCoupledSpecies< NormalForceSpecies >::ambientHumidity_, HeatFluidCoupledSpecies< NormalForceSpecies >::ambientTemperature_, HeatFluidCoupledSpecies< NormalForceSpecies >::ambientVapourConcentration_, HeatFluidCoupledSpecies< NormalForceSpecies >::evaporationCoefficientA_, HeatFluidCoupledSpecies< NormalForceSpecies >::evaporationCoefficientB_, HeatFluidCoupledSpecies< NormalForceSpecies >::latentHeatVaporization_, HeatFluidCoupledSpecies< NormalForceSpecies >::liquidDensity_, and HeatFluidCoupledSpecies< NormalForceSpecies >::massTransferCoefficient_.

◆ HeatFluidCoupledSpecies() [2/2]

◆ ~HeatFluidCoupledSpecies()

template<class NormalForceSpecies >
HeatFluidCoupledSpecies< NormalForceSpecies >::~HeatFluidCoupledSpecies
virtual

The default destructor.

199 {}

Member Function Documentation

◆ actionsAfterTimeStep()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::actionsAfterTimeStep ( BaseParticle particle) const
override
376 {
377  double dt = baseParticle->getHandler()->getDPMBase()->getTimeStep();
378  auto p = dynamic_cast<HeatFluidCoupledParticle*>(baseParticle);
379  double mass = p->getMass();
380  double surfaceArea = p->getSurfaceArea();
381  //Runge–Kutta method
382  std::array<double,2> k1 = f(p->getLiquidVolume(),p->getTemperature(), mass, surfaceArea);
383  //logger(INFO, "dL % dT %", k1[0], k1[1]);
384  std::array<double,2> k2 = f(p->getLiquidVolume()+dt*0.5*k1[0],p->getTemperature()+dt*0.5*k1[1], mass, surfaceArea);
385  std::array<double,2> k3 = f(p->getLiquidVolume()+dt*0.5*k2[0],p->getTemperature()+dt*0.5*k2[1], mass, surfaceArea);
386  std::array<double,2> k4 = f(p->getLiquidVolume()+dt*k3[0],p->getTemperature()+dt*k3[1], mass, surfaceArea);
387  double dliquidVolume = dt*(k1[0]+2*k2[0]+2*k3[0]+k4[0])/6.0;
388  double dTemperature = dt*(k1[1]+2*k2[1]+2*k3[1]+k4[1])/6.0;
389  // if liquid film volume is larger than the volume that needs to be subtracted
390  if (p->getLiquidVolume()+dliquidVolume>=0.0) {
391  p->setLiquidVolume(p->getLiquidVolume()+dliquidVolume);
392  p->setTemperature(std::max(0.0,p->getTemperature()+dTemperature));
393  //logger(INFO,"% LF % T %", p->getIndex(), p->getLiquidVolume(), p->getTemperature());
394  } else {
395  // how much to remove from liquid bridges
396  double liquidVolumeToDistribute = - (p->getLiquidVolume()+dliquidVolume);
397  p->setLiquidVolume(0.0);
398  // get liquid bridge volume
399  double liquidBridgeVolume = 0.0;
400  for (BaseInteraction* i : p-> getInteractions()) {
401  auto j = dynamic_cast<LiquidMigrationWilletInteraction*>(i);
402  liquidBridgeVolume += j->getLiquidBridgeVolume();
403  }
404  // if liquid bridge volume is larger than the volume that needs to be subtracted
405  if (liquidVolumeToDistribute<=liquidBridgeVolume) {
406  double factor = 1.0-liquidVolumeToDistribute/liquidBridgeVolume;
407  for (BaseInteraction* i : p-> getInteractions()) {
408  auto j = dynamic_cast<LiquidMigrationWilletInteraction*>(i);
409  j->setLiquidBridgeVolume(factor*j->getLiquidBridgeVolume());
410  }
411  p->setTemperature(std::max(0.0,p->getTemperature()+dTemperature));
412  } else {
413  // if both liquid bridges and liquid films are empty after the subtraction
414  if (liquidBridgeVolume!=0.0) {
415  liquidVolumeToDistribute -= liquidBridgeVolume;
416  for (BaseInteraction *i: p->getInteractions()) {
417  auto j = dynamic_cast<LiquidMigrationWilletInteraction *>(i);
418  j->setLiquidBridgeVolume(0.0);
419  }
420  }
421  double factor = 1.0+liquidVolumeToDistribute/dliquidVolume;
422  p->setTemperature(std::max(0.0,p->getTemperature()+factor*dTemperature));
423  }
424  }
425 }
Stores information about interactions between two interactable objects; often particles but could be ...
Definition: BaseInteraction.h:60
std::array< double, 2 > f(double liquidVolume, double temperature, double mass, double surfaceArea) const
f1 is used in Runge–Kutta method.
Definition: HeatFluidCoupledSpecies.h:436
Class of particles that store both temperature and liquid volume, which is adapted for the CFD-DEM st...
Definition: HeatFluidCoupledParticle.h:46
Defines the liquid bridge willet interaction between two particles or walls.
Definition: LiquidMigrationWilletInteraction.h:45
void setLiquidBridgeVolume(Mdouble liquidBridgeVolume)
Definition: LiquidMigrationWilletInteraction.cc:467
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

References BaseHandler< T >::getDPMBase(), BaseParticle::getHandler(), DPMBase::getTimeStep(), constants::i, and LiquidMigrationWilletInteraction::setLiquidBridgeVolume().

◆ f()

template<class NormalForceSpecies >
std::array< double, 2 > HeatFluidCoupledSpecies< NormalForceSpecies >::f ( double  liquidVolume_,
double  temperature_,
double  mass,
double  surfaceArea 
) const

f1 is used in Runge–Kutta method.

Computes the drying and cooling rate of a particle; based on equations in (Azmir et al., 2018), which are summarised in EvaporationModel.pdf

Parameters
liquidVolume_Liquid film volume
temperature_Temperature of particle
massMass of particle
surfaceAreaSurface area of particle
Returns
436  {
437  // Saturated vapour concentration (kg/m^3), eq(7)
438  // Dependency on temperature is valid between 1 and 200 degree Celsius
439  // (it is unphysically decreasing between 0 and 1 degree),
440  // extrapolated as linear functions
441  double dT = temperature_-273;
442  double saturatedVapourConcentration = dT < 1 ? 8.319815774e-3*temperature_/274 :
443  (((4.844e-9*dT-1.4807e-7)*dT+2.6572e-5)*dT-4.8613e-5)*dT+8.342e-3;
444 
445  // Relative activation energy of evaporation, eq(6)
446  double equilibriumActivationEnergy=-constants::R*getAmbientTemperature()*log(getAmbientHumidity());
447  // \todo what is the difference between the variable X and Y in Azmir?
448  double moistureContent = liquidVolume_*getLiquidDensity()/mass;
449  double activationEnergy=equilibriumActivationEnergy*(moistureContent-getAmbientEquilibriumMoistureContent());
450 
451  // Vapour concentrations at the particle-medium interface (kg/m^3), eq(5)
452  // Implemented such that the code does not necessarily fail if temperature is 0
453  double interfaceVapourConcentration = temperature_==0?0.0:
454  exp(-activationEnergy/(constants::R*temperature_))*saturatedVapourConcentration;
455 
456  // Drying rate of liquid film mass (kg/s), eq (4)
457  double dLiquidMass=-getMassTransferCoefficient()*surfaceArea
458  *(interfaceVapourConcentration-getAmbientVapourConcentration());
459 
460  // Drying rate of liquid film volume (kg/s), eq (3)
461  double dLiquidVolume = dLiquidMass/getLiquidDensity();
462 
463  // Heat of evaporation (J/s), eq (2)
464  double heatOfEvaporation = getLatentHeatVaporization()*(1.0+getEvaporationCoefficientA()*exp((getEvaporationCoefficientB()*getLiquidDensity()/mass)*liquidVolume_))*dLiquidMass;
465 
466  // Cooling rate of particle (K/s), eq (1)
467  double dTemperature = heatOfEvaporation / (mass * this->getHeatCapacity());
468 
469  return {dLiquidVolume, dTemperature};
470 }
Mdouble getAmbientVapourConcentration() const
Allows ambientVapourConcentration_ to be accessed.
Definition: HeatFluidCoupledSpecies.h:345
Mdouble getEvaporationCoefficientA() const
Allows evaporationCoefficientA_ to be accessed.
Definition: HeatFluidCoupledSpecies.h:284
Mdouble getLiquidDensity() const
Allows liquidDensity_ to be accessed.
Definition: HeatFluidCoupledSpecies.h:269
Mdouble getAmbientHumidity() const
Allows ambientHumidity_ to be accessed.
Definition: HeatFluidCoupledSpecies.h:314
Mdouble getAmbientEquilibriumMoistureContent() const
Allows ambientEquilibriumMoistureContent_ to be accessed.
Definition: HeatFluidCoupledSpecies.h:330
Mdouble getEvaporationCoefficientB() const
Allows evaporationCoefficientB_ to be accessed.
Definition: HeatFluidCoupledSpecies.h:299
Mdouble getAmbientTemperature() const
Allows ambientTemperature_ to be accessed.
Definition: HeatFluidCoupledSpecies.h:360
Mdouble getLatentHeatVaporization() const
Allows latentHeatVaporization_ to be accessed.
Definition: HeatFluidCoupledSpecies.h:254
Mdouble getMassTransferCoefficient() const
Allows massTransferCoefficient_ to be accessed.
Definition: HeatFluidCoupledSpecies.h:239
Mdouble getHeatCapacity() const
Allows heatCapacity_ to be accessed.
Definition: ThermalSpecies.h:128
const Mdouble R
Definition: ExtendedMath.h:50
Mdouble log(Mdouble Power)
Definition: ExtendedMath.cc:104
Mdouble exp(Mdouble Exponent)
Definition: ExtendedMath.cc:84

References mathsFunc::exp(), mathsFunc::log(), and constants::R.

◆ getAmbientEquilibriumMoistureContent()

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::getAmbientEquilibriumMoistureContent

Allows ambientEquilibriumMoistureContent_ to be accessed.

331 {
333 }

◆ getAmbientHumidity()

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::getAmbientHumidity

Allows ambientHumidity_ to be accessed.

315 {
316  return ambientHumidity_;
317 }

◆ getAmbientTemperature()

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::getAmbientTemperature

Allows ambientTemperature_ to be accessed.

361 {
362  return ambientTemperature_;
363 }

◆ getAmbientVapourConcentration()

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::getAmbientVapourConcentration

Allows ambientVapourConcentration_ to be accessed.

346 {
348 }

◆ getBaseName()

template<class NormalForceSpecies >
std::string HeatFluidCoupledSpecies< NormalForceSpecies >::getBaseName

Used in Species::getName to obtain a unique name for each Species.

234 {
235  return "HeatFluidCoupled" + NormalForceSpecies::getBaseName();
236 }

◆ getEvaporationCoefficientA()

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::getEvaporationCoefficientA

Allows evaporationCoefficientA_ to be accessed.

285 {
287 }

◆ getEvaporationCoefficientB()

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::getEvaporationCoefficientB

Allows evaporationCoefficientB_ to be accessed.

300 {
302 }

◆ getLatentHeatVaporization()

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::getLatentHeatVaporization

Allows latentHeatVaporization_ to be accessed.

255 {
257 }

◆ getLiquidDensity()

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::getLiquidDensity

Allows liquidDensity_ to be accessed.

270 {
271  return liquidDensity_;
272 }

◆ getMassTransferCoefficient()

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::getMassTransferCoefficient

Allows massTransferCoefficient_ to be accessed.

240 {
242 }

◆ read()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::read ( std::istream &  is)

Reads the species properties from an input stream.

218 {
219  std::string dummy;
221  is >> dummy >> massTransferCoefficient_;
222  is >> dummy >> latentHeatVaporization_;
223  is >> dummy >> liquidDensity_;
224  is >> dummy >> evaporationCoefficientA_;
225  is >> dummy >> evaporationCoefficientB_;
226  is >> dummy >> ambientHumidity_;
227  is >> dummy >> ambientEquilibriumMoistureContent_;
228  is >> dummy >> ambientVapourConcentration_;
229  is >> dummy >> ambientTemperature_;
230 }
void read(std::istream &is)
Reads the species properties from an input stream.
Definition: ThermalSpecies.h:112

References ThermalSpecies< NormalForceSpecies >::read().

◆ setAmbientEquilibriumMoistureContent()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::setAmbientEquilibriumMoistureContent ( Mdouble  ambientEquilibriumMoistureContent)

Allows ambientEquilibriumMoistureContent_ to be changed.

337 {
338  logger.assert_always(ambientEquilibriumMoistureContent >= 0,
339  "[HeatFluidCoupledSpecies<>::setAmbientEquilibriumMoistureContent(%)] value has to be positive",
340  ambientEquilibriumMoistureContent);
341  ambientEquilibriumMoistureContent_ = ambientEquilibriumMoistureContent;
342 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.

References logger.

◆ setAmbientHumidity()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::setAmbientHumidity ( Mdouble  ambientHumidity)

Allows ambientHumidity_ to be changed.

321 {
322  //note, this is not allowed to be zero!
323  logger.assert_always(ambientHumidity > 0,
324  "[HeatFluidCoupledSpecies<>::setAmbientHumidity(%)] value has to be positive",
325  ambientHumidity);
326  ambientHumidity_ = ambientHumidity;
327 }

References logger.

◆ setAmbientTemperature()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::setAmbientTemperature ( Mdouble  ambientTemperature)

Allows ambientTemperature_ to be changed.

367 {
368  logger.assert_always(ambientTemperature > 0,
369  "[HeatFluidCoupledSpecies<>::setAmbientTemperature(%)] value has to be positive",
370  ambientTemperature);
371  ambientTemperature_ = ambientTemperature;
372 }

References logger.

◆ setAmbientVapourConcentration()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::setAmbientVapourConcentration ( Mdouble  ambientVapourConcentration)

Allows ambientVapourConcentration_ to be changed.

352 {
353  logger.assert_always(ambientVapourConcentration >= 0,
354  "[HeatFluidCoupledSpecies<>::setAmbientVapourConcentration(%)] value has to be positive",
355  ambientVapourConcentration);
356  ambientVapourConcentration_ = ambientVapourConcentration;
357 }

References logger.

◆ setEvaporationCoefficientA()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::setEvaporationCoefficientA ( Mdouble  evaporationCoefficientA)

Allows evaporationCoefficientA_ to be changed.

291 {
292  logger.assert_always(evaporationCoefficientA >= 0,
293  "[HeatFluidCoupledSpecies<>::setEvaporationCoefficientA(%)] value has to be positive",
294  evaporationCoefficientA);
295  evaporationCoefficientA_ = evaporationCoefficientA;
296 }

References logger.

◆ setEvaporationCoefficientB()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::setEvaporationCoefficientB ( Mdouble  evaporationCoefficientB)

Allows evaporationCoefficientB_ to be changed.

306 {
307  logger.assert_always(evaporationCoefficientB <= 0,
308  "[HeatFluidCoupledSpecies<>::setEvaporationCoefficientB(%)] value has to be negative",
309  evaporationCoefficientB);
310  evaporationCoefficientB_ = evaporationCoefficientB;
311 }

References logger.

◆ setLatentHeatVaporization()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::setLatentHeatVaporization ( Mdouble  latentHeatVaporization)

Allows latentHeatVaporization_ to be changed.

261 {
262  logger.assert_always(latentHeatVaporization > 0,
263  "[HeatFluidCoupledSpecies<>::setLatentHeatVaporization(%)] value has to be positive",
264  latentHeatVaporization);
265  latentHeatVaporization_ = latentHeatVaporization;
266 }

References logger.

◆ setLiquidDensity()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::setLiquidDensity ( Mdouble  liquidDensity)

Allows liquidDensity_ to be changed.

276 {
277  logger.assert_always(liquidDensity > 0,
278  "[HeatFluidCoupledSpecies<>::setLiquidDensity(%)] value has to be positive",
279  liquidDensity);
280  liquidDensity_ = liquidDensity;
281 }

References logger.

◆ setMassTransferCoefficient()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::setMassTransferCoefficient ( Mdouble  massTransferCoefficient)

Allows massTransferCoefficient_ to be changed.

246 {
247  logger.assert_always(massTransferCoefficient > 0,
248  "[HeatFluidCoupledSpecies<>::setMassTransferCoefficient(%)] value has to be positive",
249  massTransferCoefficient);
250  massTransferCoefficient_ = massTransferCoefficient;
251 }

References logger.

◆ write()

template<class NormalForceSpecies >
void HeatFluidCoupledSpecies< NormalForceSpecies >::write ( std::ostream &  os) const

Writes the species properties to an output stream.

203 {
205  os << " massTransferCoefficient " << massTransferCoefficient_;
206  os << " latentHeatVaporization " << latentHeatVaporization_;
207  os << " liquidDensity " << liquidDensity_;
208  os << " evaporationCoefficientA " << evaporationCoefficientA_;
209  os << " evaporationCoefficientA " << evaporationCoefficientB_;
210  os << " ambientHumidity " << ambientHumidity_;
211  os << " ambientEquilibriumMoistureContent " << ambientEquilibriumMoistureContent_;
212  os << " ambientVapourConcentration " << ambientVapourConcentration_;
213  os << " ambientTemperature " << ambientTemperature_;
214 }
void write(std::ostream &os) const
Writes the species properties to an output stream.
Definition: ThermalSpecies.h:104

References ThermalSpecies< NormalForceSpecies >::write().

Member Data Documentation

◆ ambientEquilibriumMoistureContent_

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::ambientEquilibriumMoistureContent_
private

The ambient equilibrium moisture content (dimensionless, between 0 and 1).

Referenced by HeatFluidCoupledSpecies< NormalForceSpecies >::HeatFluidCoupledSpecies().

◆ ambientHumidity_

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::ambientHumidity_
private

The ambient humidity (dimensionless, between 0 and 1, but cannot be 0)

Referenced by HeatFluidCoupledSpecies< NormalForceSpecies >::HeatFluidCoupledSpecies().

◆ ambientTemperature_

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::ambientTemperature_
private

◆ ambientVapourConcentration_

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::ambientVapourConcentration_
private

The ambient vapour concentration (kg/m^3).

Referenced by HeatFluidCoupledSpecies< NormalForceSpecies >::HeatFluidCoupledSpecies().

◆ evaporationCoefficientA_

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::evaporationCoefficientA_
private

The evaporation coefficient a (dimensionless)

Referenced by HeatFluidCoupledSpecies< NormalForceSpecies >::HeatFluidCoupledSpecies().

◆ evaporationCoefficientB_

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::evaporationCoefficientB_
private

The evaporation coefficient b (dimensionless)

Referenced by HeatFluidCoupledSpecies< NormalForceSpecies >::HeatFluidCoupledSpecies().

◆ latentHeatVaporization_

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::latentHeatVaporization_
private

The latent heat of vaporization (J/kg)

Referenced by HeatFluidCoupledSpecies< NormalForceSpecies >::HeatFluidCoupledSpecies().

◆ liquidDensity_

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::liquidDensity_
private

◆ massTransferCoefficient_

template<class NormalForceSpecies >
Mdouble HeatFluidCoupledSpecies< NormalForceSpecies >::massTransferCoefficient_
private

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