CGCoordinates::XYZ Class Reference

Defines the position of the CGPoint, in the non-averaged directions, i.e. all directions on which spatial coarse-graining is applied (all directions for XYZ); all other directions are averaged over homogeneously. More...

#include <XYZ.h>

+ Inheritance diagram for CGCoordinates::XYZ:

Public Member Functions

void write (std::ostream &os) const
 Writes the coordinates in human-readable form to an ostream. More...
 
Mdouble getDistanceSquared (const Vec3D &p) const
 Returns the square of the distance between the particle p and the current CGPoint, in the non-averaged directions. More...
 
void setXYZ (Vec3D p)
 Returns the position of the current CGPoint, in the non-averaged directions. More...
 
void setHalfExtend (Vec3D halfExtend)
 Sets the extend of the current grid cell. More...
 
Mdouble getINormal (const BaseInteraction &c, const Vec3D &normal) const
 For the Interaction between particles/walls P and I, this function returns the dot product between the normal vector of the interaction and the branch vector from the current CGPoint towards I. More...
 
Mdouble getPNormal (const BaseInteraction &c, const Vec3D &normal) const
 For the Interaction between particles/walls P and I, this function returns the dot product between the normal vector of the interaction and the branch vector from the current CGPoint towards P. More...
 
Mdouble getCNormal (const BaseInteraction &c, const Vec3D &normal) const
 For the Interaction between particles/walls P and I, this function returns the dot product between the normal vector of the interaction and the branch vector from the current CGPoint towards the contact point. More...
 
Mdouble getTangentialSquared (const BaseInteraction &c, Mdouble pNormal) const
 For the Interaction between particles/walls P and I, this function returns the square of the minimum distance between the the current CGPoint and the branch vector between P and I. More...
 
- Public Member Functions inherited from CGCoordinates::BaseCoordinates
virtual Mdouble getWeight ()
 

Static Public Member Functions

static void writeNames (std::ostream &os)
 Writes the coordinate names in human-readable form to an ostream. More...
 
static Mdouble getVolumeOfAveragedDimensions (const Vec3D &min, const Vec3D &max)
 returns the factor the CGFunction has to be divided by, due to integrating the variables over the averaged dimensions, 1.0 for XYZ. More...
 
static Mdouble getLength (const Vec3D &p)
 Returns the length of the input vector in the non-averaged directions. More...
 
static Mdouble getGaussPrefactor (Mdouble width, Mdouble cutoff)
 Computes the prefactor of the Gauss CGFunction, which is dependent on the number of non-averaged dimensions. More...
 
static Mdouble getGaussIntegralPrefactor (Mdouble distance, Mdouble width, Mdouble cutoff)
 Computes the prefactor of the Gauss line integral, which is dependent on the number of non-averaged dimensions. More...
 
static void normalisePolynomialCoefficients (std::vector< Mdouble > &coefficients, Mdouble cutoff)
 Normalises the coefficients of Polynomial CGFunction such that the integral over all non-averaged dimensions is unity. More...
 
static const unsigned countVariables ()
 
static bool isResolvedIn (unsigned i)
 
static std::string getName ()
 
- Static Public Member Functions inherited from CGCoordinates::BaseCoordinates
static Mdouble getDomainVolume (const Vec3D &min, const Vec3D &max)
 

Protected Attributes

Vec3D p_
 
Vec3D pointHalfExtend_
 

Detailed Description

Defines the position of the CGPoint, in the non-averaged directions, i.e. all directions on which spatial coarse-graining is applied (all directions for XYZ); all other directions are averaged over homogeneously.

In addition to defining the spatial variable, the classes in CGCoordinates contain all functions that only depend on how many coordinate directions are locally resolved (i.e. not averaged over); e.g. getDistanceSquared and getVolumeOfAveragedDimensions.

The CGCoordinates class should be chosen such that the system in homogeneous (symmetric) in the averaged directions. E.g., periodic chutes are homogeneous in x, y and t, but varying in z, so Z should be used.

See member CGCoordinates for more details.

Member Function Documentation

◆ countVariables()

const unsigned XYZ::countVariables ( )
static

returns the number of variables (in this case three)

194 {
195  return 3;
196 }

◆ getCNormal()

Mdouble XYZ::getCNormal ( const BaseInteraction c,
const Vec3D normal 
) const

For the Interaction between particles/walls P and I, this function returns the dot product between the normal vector of the interaction and the branch vector from the current CGPoint towards the contact point.

Parameters
[in]cthe Interaction object from which iNormal is computed
Returns
cNormal, one of the three distances needed to calculate the line integral \(\psi\) which defines the stress distribution (see image).
Illustration of the line integral
118 {
119  return Vec3D::dot(c.getContactPoint() - p_, c.getNormal());
120 }
Vec3D p_
Definition: XYZ.h:164
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:56
int c
Definition: calibrate.py:100

◆ getDistanceSquared()

Mdouble XYZ::getDistanceSquared ( const Vec3D p) const

Returns the square of the distance between the particle p and the current CGPoint, in the non-averaged directions.

This function is needed to evaluate the CGFunction, as this function has the distance between the CGPoint and the Particle as an argument. To properly account for the averaging, the distance is only computed in the non-averaged directions.

Parameters
[in]pthe position of a particle for which the distance is computed.
77 {
78  return Vec3D::getLengthSquared(p_ - p);
79 }
float * p
Definition: Tutorial_Map_using.cpp:9
Mdouble getLengthSquared() const
Calculates the squared length of this Vec3D: .
Definition: Vector.cc:164

References Vec3D::getLengthSquared(), p, and p_.

◆ getGaussIntegralPrefactor()

Mdouble XYZ::getGaussIntegralPrefactor ( Mdouble  distance,
Mdouble  width,
Mdouble  cutoff 
)
static

Computes the prefactor of the Gauss line integral, which is dependent on the number of non-averaged dimensions.

The prefactor of the Gauss line integral is set such that the integral over the non-averaged dimensions is unity.

Parameters
[in]distancelength of the branch vector along which the line integral is evaluated.
[in]widthwidth (equals the standard deviation in 1D) of the Gauss CGFunction.
[in]cutoffcutoff of the Gauss CGFunction
Returns
the prefactor of the Gauss CGFunction.
161 {
162  Mdouble widthSqrt2 = width * constants::sqrt_2;
163  Mdouble a = -cutoff;
164  Mdouble b = cutoff + distance;
165  //full 2D prefactor
166  Mdouble prefactor = 1.0 / (constants::sqrt_2 * constants::sqrt_pi * width);
167  prefactor = mathsFunc::square(prefactor) / (1.0 - exp(-0.5 * mathsFunc::square(cutoff / width)));
168  return prefactor * 0.5 / (
169  +erf(b / widthSqrt2) * b
170  + widthSqrt2 / constants::sqrt_pi * exp(-mathsFunc::square(b / widthSqrt2))
171  - erf(a / widthSqrt2) * a
172  - widthSqrt2 / constants::sqrt_pi * exp(-mathsFunc::square(a / widthSqrt2))
173  );
174 }
Scalar * b
Definition: benchVecAdd.cpp:17
const Scalar * a
Definition: level2_cplx_impl.h:32
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 exp(const bfloat16 &a)
Definition: BFloat16.h:615
const Mdouble sqrt_pi
Definition: ExtendedMath.h:24
const Mdouble sqrt_2
Definition: ExtendedMath.h:26
T square(const T val)
squares a number
Definition: ExtendedMath.h:86

References a, b, Eigen::bfloat16_impl::exp(), constants::sqrt_2, constants::sqrt_pi, and mathsFunc::square().

◆ getGaussPrefactor()

Mdouble XYZ::getGaussPrefactor ( Mdouble  width,
Mdouble  cutoff 
)
static

Computes the prefactor of the Gauss CGFunction, which is dependent on the number of non-averaged dimensions.

The prefactor of the Gauss CGFunction is set such that the integral over the non-averaged dimensions is unity.

Parameters
[in]widthwidth (equals the standard deviation in 1D) of the Gauss CGFunction.
[in]cutoffcutoff of the Gauss CGFunction
Returns
the prefactor of the Gauss CGFunction.
142 {
143  //Wolfram alpha: erf(c/(sqrt(2) w))-(sqrt(2/pi) c e^(-c^2/(2 w^2)))/w
144  Mdouble prefactor = 1.0 / (constants::sqrt_2 * constants::sqrt_pi * width);
145  Mdouble cw = cutoff / width;
146  return mathsFunc::cubic(prefactor) / (
147  erf(cw / constants::sqrt_2)
149  );
150 }
T cubic(const T val)
calculates the cube of a number
Definition: ExtendedMath.h:95

References mathsFunc::cubic(), Eigen::bfloat16_impl::exp(), constants::sqrt_2, constants::sqrt_pi, and mathsFunc::square().

◆ getINormal()

Mdouble XYZ::getINormal ( const BaseInteraction c,
const Vec3D normal 
) const

For the Interaction between particles/walls P and I, this function returns the dot product between the normal vector of the interaction and the branch vector from the current CGPoint towards I.

Parameters
[in]cthe Interaction object from which iNormal is computed
Returns
iNormal, one of the three distances needed to calculate the line integral \(\psi\) which defines the stress distribution (see image).
Illustration of the line integral
98 {
99  return Vec3D::dot(c.getI()->getPosition() - p_, c.getNormal());
100 }

References calibrate::c, Vec3D::dot(), and p_.

◆ getLength()

Mdouble XYZ::getLength ( const Vec3D p)
static

Returns the length of the input vector in the non-averaged directions.

Parameters
[in]pvector whose length should be determined
Returns
length of the vector in the non-averaged directions
Todo:
87 {
88  return sqrt(p.X * p.X + p.Y * p.Y + p.Z * p.Z);
89 }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134

References p, and sqrt().

◆ getName()

std::string XYZ::getName ( )
static
199 {
200  return "XYZ";
201 }

◆ getPNormal()

Mdouble XYZ::getPNormal ( const BaseInteraction c,
const Vec3D normal 
) const

For the Interaction between particles/walls P and I, this function returns the dot product between the normal vector of the interaction and the branch vector from the current CGPoint towards P.

Parameters
[in]cthe Interaction object from which iNormal is computed
Returns
pNormal, one of the three distances needed to calculate the line integral \(\psi\) which defines the stress distribution (see image).
Illustration of the line integral
108 {
109  return Vec3D::dot(c.getP()->getPosition() - p_, c.getNormal());
110 }

◆ getTangentialSquared()

Mdouble XYZ::getTangentialSquared ( const BaseInteraction c,
Mdouble  pNormal 
) const

For the Interaction between particles/walls P and I, this function returns the square of the minimum distance between the the current CGPoint and the branch vector between P and I.

Parameters
[in]cthe Interaction object from which iNormal is computed
[in]pNormalthe output of getPNormal needed for the computation.
Returns
iNormal, one of the three distances needed to calculate the line integral \(\psi\) which defines the stress distribution (see image).
Illustration of the line integral
130 {
131  return Vec3D::getLengthSquared(c.getP()->getPosition() - p_) - mathsFunc::square(pNormal);
132 }

References calibrate::c, Vec3D::getLengthSquared(), p_, and mathsFunc::square().

◆ getVolumeOfAveragedDimensions()

Mdouble XYZ::getVolumeOfAveragedDimensions ( const Vec3D min,
const Vec3D max 
)
static

returns the factor the CGFunction has to be divided by, due to integrating the variables over the averaged dimensions, 1.0 for XYZ.

If averaged dimensions are present (i.e. for all Coordinates except XYZ), the CGFunction has to be divided by a factor (here called volume) due to integrating the variables over the averaged dimensions.

Parameters
[in]minthe lower limits of the mesh domain (xMin, yMin, zMin)
[in]maxthe upper limits of the mesh domain (xMax, yMax, zMax)
Returns
the volume factor
57 {
58  return 1.0;
59 }

◆ isResolvedIn()

static bool CGCoordinates::XYZ::isResolvedIn ( unsigned  i)
inlinestatic
155 {return true;}

◆ normalisePolynomialCoefficients()

void XYZ::normalisePolynomialCoefficients ( std::vector< Mdouble > &  coefficients,
Mdouble  cutoff 
)
static

Normalises the coefficients of Polynomial CGFunction such that the integral over all non-averaged dimensions is unity.

The volume is computed as

\[volume= \int_0^1\sum_{i=1}^n c_i r/c^i 4 pi r^2 dr = 4 pi \sum_{i=1}^n c_i/(i+3) \]

with 4 pi r^2 the surface area of a sphere.

Parameters
[in,out]coefficientsthe coefficients of Polynomial CGFunctions.
[in]cutoffcutoff of the Gauss CGFunction
184 {
185  Mdouble volume = 0.0;
186  for (std::size_t i = 0; i < coefficients.size(); i++)
187  volume += coefficients[i] / static_cast<Mdouble>(i + 3);
188  volume *= 4.0 * constants::pi * mathsFunc::cubic(cutoff);
189  for (double& coefficient : coefficients)
190  coefficient /= volume;
191 }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const Mdouble pi
Definition: ExtendedMath.h:23

References mathsFunc::cubic(), i, and constants::pi.

◆ setHalfExtend()

void CGCoordinates::XYZ::setHalfExtend ( Vec3D  halfExtend)
inline

Sets the extend of the current grid cell.

102 { pointHalfExtend_ = halfExtend; }
Vec3D pointHalfExtend_
Definition: XYZ.h:165

References pointHalfExtend_.

◆ setXYZ()

void XYZ::setXYZ ( Vec3D  p)

Returns the position of the current CGPoint, in the non-averaged directions.

Parameters
[in]pthe position that is to be set.
65 {
66  p_ = p;
67 }

References p, and p_.

Referenced by VolumeCoupling::getProjectionMatrix().

◆ write()

void XYZ::write ( std::ostream &  os) const

Writes the coordinates in human-readable form to an ostream.

Parameters
[out]osthe ostream file to which the position is written
20 {
21  os << p_ << ' ';
22 }

References p_.

◆ writeNames()

void XYZ::writeNames ( std::ostream &  os)
static

Writes the coordinate names in human-readable form to an ostream.

12 {
13  os << "x y z ";
14 }

Member Data Documentation

◆ p_

Vec3D CGCoordinates::XYZ::p_
protected

The position of the current CGPoint.

Referenced by getDistanceSquared(), getINormal(), getTangentialSquared(), setXYZ(), and write().

◆ pointHalfExtend_

Vec3D CGCoordinates::XYZ::pointHalfExtend_
protected

Referenced by setHalfExtend().


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