CGFields::StandardFields Class Reference

Contains the computed field values, like density, momentum and stress. More...

#include <StandardFields.h>

+ Inheritance diagram for CGFields::StandardFields:

Public Member Functions

 StandardFields ()
 Default constructor, sets all field values to zero. More...
 
 StandardFields (const StandardFields &P)=default
 Default copy constructor, copies the values of all fields. More...
 
 ~StandardFields ()=default
 Destructor, it simply destructs the StandardFields and all the objects it contains. More...
 
void write (std::ostream &os) const
 Writes class content into an output stream, typically a stat file. More...
 
void output (std::ostream &os) const
 Writes human-readable class content into an output stream, typically a stat file. More...
 
void setZero ()
 Sets all fields to zero. More...
 
StandardFields getSquared () const
 Returns the square of all field values (to calculate standard deviation). More...
 
StandardFieldsoperator= (const StandardFields &P)
 Copies all field values. More...
 
StandardFieldsoperator+= (const StandardFields &P)
 Adds the field values on the RHS to the LHS of the equation. More...
 
StandardFieldsoperator-= (const StandardFields &P)
 Subtracts the field values on the RHS from the LHS of the equation. More...
 
StandardFieldsoperator/= (Mdouble a)
 Divides the field values on the LHS by the RHS of the equation. More...
 
StandardFields operator* (Mdouble a) const
 Multiplies the field values on the left of the '*' by the scalar value on the right of the '*' and returns the answer. More...
 
void addParticleStatistics (Mdouble phi, const StandardFields &currentInteraction)
 This function should be called from within a loop over all particles to compute all the fields that are defined as a sum over all particles (e.g. density, momentum). More...
 
void addParticleDifferentialStatistics (Vec3D &dphi, const StandardFields &currentInteraction)
 
void addInteractionStatistics (Mdouble psi, const StandardFields &currentInteraction)
 This function should be called from within a loop over all Interactions to compute all the fields that are defined as a sum over all Interactions (e.g. stress). More...
 
void addContactPointStatistics (Mdouble phi, const StandardFields &currentInteraction)
 This function should be called from within a loop over all Interactions to compute all the fields that are defined as a sum over all Interactions with external objects (e.g. IFD). More...
 
void setFields (const BaseInteraction &c, IntegralType type)
 
void setCylindricalFields (const BaseInteraction &c, IntegralType type)
 
void setFields (const BaseParticle &p)
 
void setCylindricalFields (const BaseParticle &p)
 
Mdouble getVolumeFraction () const
 
Mdouble getDensity () const
 
Vec3D getMomentum () const
 
MatrixSymmetric3D getMomentumFlux () const
 
Matrix3D getContactStress () const
 
Vec3D getInteractionForceDensity () const
 
Mdouble getParticleSizeDensity (size_t i) const
 
std::array< Mdouble, 6 > getParticleSizeDensity () const
 
std::array< Mdouble, 6 > getParticleSizeMomenta () const
 
std::array< Mdouble, 6 > getCentralParticleSizeMomenta () const
 
std::array< Mdouble, 6 > getStandardisedParticleSizeMomenta () const
 
void outputStandardisedParticleSizeMomenta (std::ostream &os) const
 

Static Public Member Functions

static void writeNames (std::ostream &os, unsigned countVariables)
 
static bool doInteractionStatistics ()
 Returns true if the class contains fields that are defined as a sum over all Interactions (e.g. stress), else returns false. More...
 
static bool evaluateFixedParticles ()
 
static bool isDifferentialField ()
 

Private Attributes

Mdouble volumeFraction_
 
Mdouble density_
 
Vec3D momentum_
 
MatrixSymmetric3D momentumFlux_
 
Matrix3D contactStress_
 
Vec3D interactionForceDensity_
 
std::array< Mdouble, 6 > particleSizeDensity_
 

Detailed Description

Contains the computed field values, like density, momentum and stress.

CGPoints inherits from this class; CGPoints::evaluate adds to the values of these variables.

Todo:
These are currently the only fields that are computed. However, this class is destined to be extended to contain additional information such as fabric, energy, local angular momentum. Also, a simpler version is planned, where only particle statistics are evaluated (density and momentum).

Constructor & Destructor Documentation

◆ StandardFields() [1/2]

CGFields::StandardFields::StandardFields ( )

Default constructor, sets all field values to zero.

33 {
34  setZero();
35 #ifdef DEBUG_CONSTRUCTOR
36  std::cerr << "StandardFields::StandardFields() finished" << std::endl;
37 #endif
38 }
void setZero()
Sets all fields to zero.
Definition: StandardFields.cc:84

References setZero().

◆ StandardFields() [2/2]

CGFields::StandardFields::StandardFields ( const StandardFields P)
default

Default copy constructor, copies the values of all fields.

◆ ~StandardFields()

CGFields::StandardFields::~StandardFields ( )
default

Destructor, it simply destructs the StandardFields and all the objects it contains.

Member Function Documentation

◆ addContactPointStatistics()

void CGFields::StandardFields::addContactPointStatistics ( Mdouble  phi,
const StandardFields currentInteraction 
)

This function should be called from within a loop over all Interactions to compute all the fields that are defined as a sum over all Interactions with external objects (e.g. IFD).

Parameters
[in]phithe value of the cg function for the contact point of c and the current CGPoint
[in]cthe interaction which is used in the cg function
218 {
219  interactionForceDensity_ += currentInteraction.getInteractionForceDensity() * phi;
220 }
Vec3D interactionForceDensity_
Definition: StandardFields.h:280

References getInteractionForceDensity(), and interactionForceDensity_.

◆ addInteractionStatistics()

void CGFields::StandardFields::addInteractionStatistics ( Mdouble  psi,
const StandardFields currentInteraction 
)

This function should be called from within a loop over all Interactions to compute all the fields that are defined as a sum over all Interactions (e.g. stress).

Parameters
[in]psithe value of the line integral from C to P at the current CGPoint
[in]cthe contact which is used in the line integral
208 {
209  contactStress_ += currentInteraction.getContactStress() * psi;
210 }
Matrix3D contactStress_
Definition: StandardFields.h:270

References contactStress_, and getContactStress().

◆ addParticleDifferentialStatistics()

void CGFields::StandardFields::addParticleDifferentialStatistics ( Vec3D dphi,
const StandardFields currentInteraction 
)
200 {
201 }

◆ addParticleStatistics()

void CGFields::StandardFields::addParticleStatistics ( Mdouble  phi,
const StandardFields currentInteraction 
)

This function should be called from within a loop over all particles to compute all the fields that are defined as a sum over all particles (e.g. density, momentum).

Parameters
[in]phithe value of the cg function at the current CGPoint
[in]currentInteractionthe fields which are produced due to the particle at its centre
190 {
191  volumeFraction_ += currentInteraction.getVolumeFraction() * phi;
192  density_ += currentInteraction.getDensity() * phi;
193  momentum_ += currentInteraction.getMomentum() * phi;
194  momentumFlux_ += currentInteraction.getMomentumFlux() * phi;
195  for (size_t i = 0; i < particleSizeDensity_.size(); ++i)
196  particleSizeDensity_[i] += currentInteraction.getParticleSizeDensity(i) * phi;
197 }
Vec3D momentum_
Definition: StandardFields.h:246
MatrixSymmetric3D momentumFlux_
Definition: StandardFields.h:258
Mdouble density_
Definition: StandardFields.h:234
std::array< Mdouble, 6 > particleSizeDensity_
Definition: StandardFields.h:291
Mdouble volumeFraction_
Definition: StandardFields.h:226
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

References density_, getDensity(), getMomentum(), getMomentumFlux(), getParticleSizeDensity(), getVolumeFraction(), constants::i, momentum_, momentumFlux_, particleSizeDensity_, and volumeFraction_.

◆ doInteractionStatistics()

bool CGFields::StandardFields::doInteractionStatistics ( )
static

Returns true if the class contains fields that are defined as a sum over all Interactions (e.g. stress), else returns false.

If the functions returns false, addInteractionStatistics and addContactPointStatistics are never called.

Returns
True if the class contains fields that are defined as a sum over all Interactions, else false.
229 {
230  return true;
231 }

◆ evaluateFixedParticles()

static bool CGFields::StandardFields::evaluateFixedParticles ( )
inlinestatic
205  {
206  return false;
207  }

◆ getCentralParticleSizeMomenta()

std::array< Mdouble, 6 > CGFields::StandardFields::getCentralParticleSizeMomenta ( ) const
303 {
304  //https://en.wikipedia.org/wiki/Central_moment
305  auto momenta = getParticleSizeMomenta();
306  Mdouble mean = momenta[1];
307  momenta[5] += -5 * mean * momenta[4] + 10 * mean * mean * momenta[3]
308  - 10 * mean * mean * mean * momenta[2] + 4 * mean * mean * mean * mean * mean;
309  momenta[4] += -4 * mean * momenta[3] + 6 * mean * mean * momenta[2] - 3 * mean * mean * mean * mean;
310  momenta[3] += -3 * mean * momenta[2] + 2 * mean * mean * mean;
311  momenta[2] += -mean * mean;
312  return momenta;
313 }
double Mdouble
Definition: GeneralDefine.h:34
std::array< Mdouble, 6 > getParticleSizeMomenta() const
Definition: StandardFields.cc:291

References getParticleSizeMomenta().

Referenced by getStandardisedParticleSizeMomenta().

◆ getContactStress()

Matrix3D CGFields::StandardFields::getContactStress ( ) const
inline
177  {
178  return contactStress_;
179  }

References contactStress_.

Referenced by addInteractionStatistics().

◆ getDensity()

Mdouble CGFields::StandardFields::getDensity ( ) const
inline
162  {
163  return density_;
164  }

References density_.

Referenced by addParticleStatistics(), and testCGHandler().

◆ getInteractionForceDensity()

Vec3D CGFields::StandardFields::getInteractionForceDensity ( ) const
inline
182  {
184  }

References interactionForceDensity_.

Referenced by addContactPointStatistics(), and testCGHandler().

◆ getMomentum()

Vec3D CGFields::StandardFields::getMomentum ( ) const
inline
167  {
168  return momentum_;
169  }

References momentum_.

Referenced by addParticleStatistics().

◆ getMomentumFlux()

MatrixSymmetric3D CGFields::StandardFields::getMomentumFlux ( ) const
inline
172  {
173  return momentumFlux_;
174  }

References momentumFlux_.

Referenced by addParticleStatistics().

◆ getParticleSizeDensity() [1/2]

std::array<Mdouble, 6> CGFields::StandardFields::getParticleSizeDensity ( ) const
inline
192  {
193  return particleSizeDensity_;
194  }

References particleSizeDensity_.

◆ getParticleSizeDensity() [2/2]

Mdouble CGFields::StandardFields::getParticleSizeDensity ( size_t  i) const
inline
187  {
188  return particleSizeDensity_[i];
189  }

References constants::i, and particleSizeDensity_.

Referenced by addParticleStatistics().

◆ getParticleSizeMomenta()

std::array< Mdouble, 6 > CGFields::StandardFields::getParticleSizeMomenta ( ) const
292 {
293  //https://en.wikipedia.org/wiki/Moment_(mathematics)
294  auto momenta = particleSizeDensity_;
295  for (size_t i = 1; i < particleSizeDensity_.size(); ++i)
296  {
297  momenta[i] /= particleSizeDensity_[0];
298  }
299  return momenta;
300 }

References constants::i, and particleSizeDensity_.

Referenced by getCentralParticleSizeMomenta().

◆ getSquared()

StandardFields CGFields::StandardFields::getSquared ( ) const

Returns the square of all field values (to calculate standard deviation).

Returns
a CGField containing the square of the values in the current object
99 {
101  P.volumeFraction_ = mathsFunc::square(volumeFraction_);
102  P.density_ = mathsFunc::square(density_);
103  P.momentum_ = Vec3D::square(momentum_);
104  P.momentumFlux_ = MatrixSymmetric3D::square(momentumFlux_);
105  P.contactStress_ = Matrix3D::square(contactStress_);
106  P.interactionForceDensity_ = Vec3D::square(interactionForceDensity_);
107  P.particleSizeDensity_ = particleSizeDensity_;
108  for (auto& p : P.particleSizeDensity_) p *= p;
109  return P;
110 }
StandardFields()
Default constructor, sets all field values to zero.
Definition: StandardFields.cc:32
static Matrix3D square(const Matrix3D &A)
Calculates the pointwise square.
Definition: Matrix.cc:298
static MatrixSymmetric3D square(const MatrixSymmetric3D &A)
Calculates the pointwise square.
Definition: MatrixSymmetric.cc:238
static Vec3D square(const Vec3D &a)
Calculates the pointwise square of a Vec3D.
Definition: Vector.cc:114
double P
Uniform pressure.
Definition: TwenteMeshGluing.cpp:73
T square(const T val)
squares a number
Definition: ExtendedMath.h:106

References contactStress_, density_, interactionForceDensity_, momentum_, momentumFlux_, Global_Physical_Variables::P, particleSizeDensity_, Matrix3D::square(), MatrixSymmetric3D::square(), mathsFunc::square(), Vec3D::square(), and volumeFraction_.

◆ getStandardisedParticleSizeMomenta()

std::array< Mdouble, 6 > CGFields::StandardFields::getStandardisedParticleSizeMomenta ( ) const
316 {
317  auto momenta = getCentralParticleSizeMomenta();
318  Mdouble std = std::sqrt(momenta[2]);
319  momenta[3] /= std * std * std;
320  momenta[4] /= std * std * std * std;
321  momenta[5] /= std * std * std * std * std;
322  return momenta;
323 }
std::array< Mdouble, 6 > getCentralParticleSizeMomenta() const
Definition: StandardFields.cc:302

References getCentralParticleSizeMomenta().

Referenced by outputStandardisedParticleSizeMomenta().

◆ getVolumeFraction()

Mdouble CGFields::StandardFields::getVolumeFraction ( ) const
inline
157  {
158  return volumeFraction_;
159  }

References volumeFraction_.

Referenced by addParticleStatistics().

◆ isDifferentialField()

static bool CGFields::StandardFields::isDifferentialField ( )
inlinestatic

A bool that determines if the derivative of the CG function has to be computed

214  {
215  return false;
216  }

◆ operator*()

StandardFields CGFields::StandardFields::operator* ( Mdouble  a) const

Multiplies the field values on the left of the '*' by the scalar value on the right of the '*' and returns the answer.

Parameters
[in]athe scalar that we multiply with
Returns
the CGField to which the multiplied values are written
156 {
157  StandardFields p;
158  p.volumeFraction_ = volumeFraction_ * a;
159  p.density_ = density_ * a;
160  p.momentum_ = momentum_ * a;
161  p.momentumFlux_ = momentumFlux_ * a;
162  p.contactStress_ = contactStress_ * a;
163  p.interactionForceDensity_ = interactionForceDensity_ * a;
164  for (size_t i = 0; i < particleSizeDensity_.size(); ++i)
165  p.particleSizeDensity_[i] = particleSizeDensity_[i] * a;
166  return p;
167 }

References contactStress_, density_, constants::i, interactionForceDensity_, momentum_, momentumFlux_, particleSizeDensity_, and volumeFraction_.

◆ operator+=()

StandardFields & CGFields::StandardFields::operator+= ( const StandardFields P)

Adds the field values on the RHS to the LHS of the equation.

Parameters
[in]Pthe CGField that has to be added
Returns
the CGField to which the values are added
124 {
125  volumeFraction_ += P.volumeFraction_;
126  density_ += P.density_;
127  momentum_ += P.momentum_;
128  momentumFlux_ += P.momentumFlux_;
129  contactStress_ += P.contactStress_;
130  interactionForceDensity_ += P.interactionForceDensity_;
131  for (size_t i = 0; i < particleSizeDensity_.size(); ++i) particleSizeDensity_[i] += P.particleSizeDensity_[i];
132  return *this;
133 }

References contactStress_, density_, constants::i, interactionForceDensity_, momentum_, momentumFlux_, Global_Physical_Variables::P, particleSizeDensity_, and volumeFraction_.

◆ operator-=()

StandardFields & CGFields::StandardFields::operator-= ( const StandardFields P)

Subtracts the field values on the RHS from the LHS of the equation.

Parameters
[in]Pthe CGField that has to be subtracted
Returns
the CGField from which the values are subtracted
140 {
141  volumeFraction_ -= P.volumeFraction_;
142  density_ -= P.density_;
143  momentum_ -= P.momentum_;
144  momentumFlux_ -= P.momentumFlux_;
145  contactStress_ -= P.contactStress_;
146  interactionForceDensity_ -= P.interactionForceDensity_;
147  for (size_t i = 0; i < particleSizeDensity_.size(); ++i) particleSizeDensity_[i] -= P.particleSizeDensity_[i];
148  return *this;
149 }

References contactStress_, density_, constants::i, interactionForceDensity_, momentum_, momentumFlux_, Global_Physical_Variables::P, particleSizeDensity_, and volumeFraction_.

◆ operator/=()

StandardFields & CGFields::StandardFields::operator/= ( Mdouble  a)

Divides the field values on the LHS by the RHS of the equation.

Parameters
[in]athe scalar that we divide by
Returns
the CGField to which the divided values are written
174 {
175  volumeFraction_ /= a;
176  density_ /= a;
177  momentum_ /= a;
178  momentumFlux_ /= a;
179  contactStress_ /= a;
181  for (auto& p : particleSizeDensity_) p /= a;
182  return *this;
183 }

References contactStress_, density_, interactionForceDensity_, momentum_, momentumFlux_, particleSizeDensity_, and volumeFraction_.

◆ operator=()

StandardFields & CGFields::StandardFields::operator= ( const StandardFields P)
default

Copies all field values.

Parameters
[in]Pthe CGField that has to be copied
Returns
the CGField into which the values are copied

◆ output()

void CGFields::StandardFields::output ( std::ostream &  os) const

Writes human-readable class content into an output stream, typically a stat file.

Parameters
[in,out]osthe ostream into which the data is written.
73 {
74  os << "VolumeFraction " << volumeFraction_;
75  os << " Density " << density_;
76  os << " Momentum " << momentum_;
77  os << " MomentumFlux " << momentumFlux_;
78  os << " ContactStress " << contactStress_;
79  os << " InteractionForce " << interactionForceDensity_;
80  os << " ParticleSize";
81  for (const auto p : particleSizeDensity_) os << ' ' << p;
82 }

References contactStress_, density_, interactionForceDensity_, momentum_, momentumFlux_, particleSizeDensity_, and volumeFraction_.

◆ outputStandardisedParticleSizeMomenta()

void CGFields::StandardFields::outputStandardisedParticleSizeMomenta ( std::ostream &  os) const
280 {
281  //https://en.wikipedia.org/wiki/Moment_%28mathematics%29#Higher_moments
282  auto momenta = getStandardisedParticleSizeMomenta();
283  os << "Particle number: " << momenta[0] << '\n';
284  os << "Mean radius: " << momenta[1] << '\n';
285  os << "Standard deviation: " << std::sqrt(momenta[2]) << '\n';
286  os << "Skewness: " << momenta[3] << '\n'; //left- or right-sided
287  os << "Kurtosis: " << momenta[4] << '\n'; //heavy or light tailed (3 for normal distribution)
288  os << "5-th moment: " << momenta[5] << '\n'; //hyper skewness (0 for normal dist)
289 }
std::array< Mdouble, 6 > getStandardisedParticleSizeMomenta() const
Definition: StandardFields.cc:315

References getStandardisedParticleSizeMomenta().

◆ setCylindricalFields() [1/2]

void CGFields::StandardFields::setCylindricalFields ( const BaseInteraction c,
IntegralType  type 
)
266 {
267  setFields(c, type);
270 }
const Vec3D & getContactPoint() const
Gets constant reference to contact point (vector).
Definition: BaseInteraction.h:234
void setFields(const BaseInteraction &c, IntegralType type)
Definition: StandardFields.cc:233
Matrix3D getCylindricalTensorField(const Vec3D &p) const
Returns the matrix in cylindrical coordinates.
Definition: Matrix.cc:394
Vec3D getCylindricalTensorField(const Vec3D &position) const
Returns this vector field at point p to cylindrical coordinates.
Definition: Vector.cc:261

References contactStress_, BaseInteraction::getContactPoint(), Matrix3D::getCylindricalTensorField(), Vec3D::getCylindricalTensorField(), interactionForceDensity_, and setFields().

◆ setCylindricalFields() [2/2]

void CGFields::StandardFields::setCylindricalFields ( const BaseParticle p)
273 {
274  setFields(p);
277 }
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218
MatrixSymmetric3D getCylindricalTensorField(const Vec3D &p) const
Returns the matrix in cylindrical coordinates.
Definition: MatrixSymmetric.cc:326

References MatrixSymmetric3D::getCylindricalTensorField(), Vec3D::getCylindricalTensorField(), BaseInteractable::getPosition(), momentum_, momentumFlux_, and setFields().

◆ setFields() [1/2]

void CGFields::StandardFields::setFields ( const BaseInteraction c,
IntegralType  type 
)
234 {
235  //P is always real, but I not
236  if (type == IntegralType::I_TO_P)
237  {
239  }
240  else if (type == IntegralType::I_TO_CONTACT)
241  {
243  }
244  else
245  {
247  }
249 }
const Vec3D & getForce() const
Gets the current force (vector) between the two interacting objects.
Definition: BaseInteraction.h:210
Vec3D getIC() const
Definition: BaseInteraction.cc:466
Vec3D getCP() const
Definition: BaseInteraction.cc:472
Vec3D getIP() const
Definition: BaseInteraction.cc:459
static Matrix3D dyadic(const Vec3D &a, const Vec3D &b)
Calculates the dyadic product of a two Vec3D: .
Definition: Matrix.cc:323

References contactStress_, Matrix3D::dyadic(), BaseInteraction::getCP(), BaseInteraction::getForce(), BaseInteraction::getIC(), BaseInteraction::getIP(), I_TO_CONTACT, I_TO_P, and interactionForceDensity_.

Referenced by setCylindricalFields().

◆ setFields() [2/2]

void CGFields::StandardFields::setFields ( const BaseParticle p)
252 {
254  density_ = p.getMass();
255  momentum_ = p.getVelocity() * p.getMass();
257  Mdouble particleSize = 1.0;
258  for (auto& ps : particleSizeDensity_)
259  {
260  ps = particleSize;
261  particleSize *= p.getRadius();
262  }
263 }
virtual const Vec3D & getVelocity() const
Returns the velocity of this interactable.
Definition: BaseInteractable.cc:329
Mdouble getRadius() const
Returns the particle's radius.
Definition: BaseParticle.h:348
virtual Mdouble getVolume() const
Get Particle volume function, which required a reference to the Species vector. It returns the volume...
Definition: BaseParticle.cc:143
Mdouble getMass() const
Returns the particle's mass.
Definition: BaseParticle.h:322
static MatrixSymmetric3D selfDyadic(const Vec3D &a)
Calculates the dyadic product of a Vec3D with itself: .
Definition: MatrixSymmetric.cc:263

References density_, BaseParticle::getMass(), BaseParticle::getRadius(), BaseInteractable::getVelocity(), BaseParticle::getVolume(), momentum_, momentumFlux_, particleSizeDensity_, MatrixSymmetric3D::selfDyadic(), and volumeFraction_.

◆ setZero()

void CGFields::StandardFields::setZero ( )

Sets all fields to zero.

85 {
86  volumeFraction_ = 0.0;
87  density_ = 0.0;
92  for (auto& p : particleSizeDensity_) p = 0;
93 }
void setZero()
Sets all elements to zero.
Definition: Matrix.cc:75
void setZero()
Sets all elements to zero.
Definition: MatrixSymmetric.cc:70
void setZero()
Sets all elements to zero.
Definition: Vector.cc:43

References contactStress_, density_, interactionForceDensity_, momentum_, momentumFlux_, particleSizeDensity_, Matrix3D::setZero(), MatrixSymmetric3D::setZero(), Vec3D::setZero(), and volumeFraction_.

Referenced by StandardFields().

◆ write()

void CGFields::StandardFields::write ( std::ostream &  os) const

Writes class content into an output stream, typically a stat file.

Parameters
[in,out]osthe ostream into which the data is written.
59 {
60  os << volumeFraction_;
61  os << ' ' << density_;
62  os << ' ' << momentum_;
63  os << ' ' << momentumFlux_;
64  os << ' ' << contactStress_;
65  os << ' ' << interactionForceDensity_;
66  for (const auto p : particleSizeDensity_) os << ' ' << p;
67 }

References contactStress_, density_, interactionForceDensity_, momentum_, momentumFlux_, particleSizeDensity_, and volumeFraction_.

◆ writeNames()

void CGFields::StandardFields::writeNames ( std::ostream &  os,
unsigned  countVariables 
)
static
Parameters
[in,out]osthe ostream into which the data is written.
[in]countVariablesThe number of variables in the field (including time), e.g. 1 for O, 4 for XYZ
45 {
46  os << countVariables + 1 << ":VolumeFraction "; //volume
47  os << countVariables + 2 << ":Density "; //mass
48  os << countVariables + 3 << '-' << countVariables + 5 << ":Momentum ";
49  os << countVariables + 6 << '-' << countVariables + 11 << ":MomentumFlux ";
50  os << countVariables + 12 << '-' << countVariables + 20 << ":ContactStress ";
51  os << countVariables + 21 << '-' << countVariables + 23 << ":InteractionForce ";
52  os << countVariables + 24 << '-' << countVariables + 29 << ":ParticleSize ";
53 }

Member Data Documentation

◆ contactStress_

Matrix3D CGFields::StandardFields::contactStress_
private

Contact stress, computed as the sum over all contacts between particle i and particle/wall/fixed particle j

\[\mathbf \sigma^{c}(\vec r,t)=\sum_{ij} \vec f_{ij} \otimes \vec l_{ij} \psi(\vec r,\vec r_i,\vec r_j),\]

with contact force \(\vec f_{ij}\), branch vector \(\vec l_{ij}= \vec c_{ij}-\vec r_i\), particle position \(\vec r_i\), contact point \(\vec c_{ij}\), and cg line integral \(\psi(\vec r,\vec r_i,\vec r_j)\), see CGFunctions::Gauss::evaluateCGIntegral.

Referenced by addInteractionStatistics(), getContactStress(), getSquared(), operator*(), operator+=(), operator-=(), operator/=(), output(), setCylindricalFields(), setFields(), setZero(), and write().

◆ density_

Mdouble CGFields::StandardFields::density_
private

(Mass) density, computed as the sum over all particles i

\[\rho(\vec r,t)=\sum_i m_i \phi(\vec r,\vec r_i)\]

with particle mass m_i and cg function \(\phi(\vec r,\vec r_i),\), see CGFunctions::Gauss::evaluateCGFunction.

Referenced by addParticleStatistics(), getDensity(), getSquared(), operator*(), operator+=(), operator-=(), operator/=(), output(), setFields(), setZero(), and write().

◆ interactionForceDensity_

Vec3D CGFields::StandardFields::interactionForceDensity_
private

Interaction force density, computed as the sum over all contacts between particle i and (external) wall/fixed particle j

\[\vec{IFD}(\vec r,t)=\sum_{ij} f_{ij} \phi(\vec r,\vec c_{ij}),\]

with contact force \(\vec f_{ij}\) and cg function \(\phi(\vec r,\vec r_i)\), see CGFunctions::Gauss::evaluateCGFunction.

Referenced by addContactPointStatistics(), getInteractionForceDensity(), getSquared(), operator*(), operator+=(), operator-=(), operator/=(), output(), setCylindricalFields(), setFields(), setZero(), and write().

◆ momentum_

Vec3D CGFields::StandardFields::momentum_
private

Momentum, computed as the sum over all particles i

\[\vec j(\vec r,t)=\sum_i m_i \vec v_i\phi(\vec r,\vec r_i),\]

with particle momentum \(m_i\vec v_i\) and cg function \(\phi(\vec r,\vec r_i)\), see CGFunctions::Gauss::evaluateCGFunction.

Velocity can be calculated in a post-processing step from momentum and density as

\[\vec V(\vec r,t) = \frac{\vec j(\vec r,t)}{\rho(\vec r,t)}.\]

Referenced by addParticleStatistics(), getMomentum(), getSquared(), operator*(), operator+=(), operator-=(), operator/=(), output(), setCylindricalFields(), setFields(), setZero(), and write().

◆ momentumFlux_

MatrixSymmetric3D CGFields::StandardFields::momentumFlux_
private

Momentum, computed as the sum over all particles i

\[\mathbf k(\vec r,t)=\sum_i m_i \vec v_i \otimes \vec v_i\phi(\vec r,\vec r_i),\]

with particle momentum flux \(m_i\vec v_i \otimes \vec v_i\) and cg function \(\phi(\vec r,\vec r_i)\), see CGFunctions::Gauss::evaluateCGFunction.

Kinetic stress can be calculated in a post-processing step from momentum flux and density as

\[\mathbf \sigma^{k} = \mathbf k - \rho \vec V \otimes \vec V.\]

Referenced by addParticleStatistics(), getMomentumFlux(), getSquared(), operator*(), operator+=(), operator-=(), operator/=(), output(), setCylindricalFields(), setFields(), setZero(), and write().

◆ particleSizeDensity_

std::array<Mdouble, 6> CGFields::StandardFields::particleSizeDensity_
private

Density of particle size, and powers of particle size,

\[ \vec{PS}_k(\vec r,t)= \sum_i r_i^k \phi(\vec r),\]

with radius \(r_i\).

Used to compute the first five moments of the particle size distribution,

\[ <r^k>(\vec r,t)= \vec{PS}_k / \vec{PS}_0, \]

with \(\vec{PS}_0\) denoting the number density.

Referenced by addParticleStatistics(), getParticleSizeDensity(), getParticleSizeMomenta(), getSquared(), operator*(), operator+=(), operator-=(), operator/=(), output(), setFields(), setZero(), and write().

◆ volumeFraction_

Mdouble CGFields::StandardFields::volumeFraction_
private

Particle volume fraction, computed as the sum over all particles i

\[\nu(\vec r,t)=\sum_i V_i \phi(\vec r,\vec r_i),\]

with particle volume V_i and cg function \(\phi(\vec r,\vec r_i)\), see CGFunctions::Gauss::evaluateCGFunction.

Referenced by addParticleStatistics(), getSquared(), getVolumeFraction(), operator*(), operator+=(), operator-=(), operator/=(), output(), setFields(), setZero(), and write().


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