Vec3D Class Reference

#include <Vector.h>

Public Member Functions

 Vec3D ()
 constructor More...
 
 Vec3D (const SmallVector< 3 > &vector)
 
 Vec3D (const Mdouble x, const Mdouble y, const Mdouble z)
 Alternative constructor, taking the three elements as arguments. More...
 
 Vec3D (std::array< double, 3 > a)
 
void setZero ()
 Sets all elements to zero. More...
 
void setNaN ()
 Sets all elements to NaN. More...
 
bool isZero () const
 Checks if ALL elements are zero. More...
 
bool isNaN () const
 Checks if ALL elements are zero. More...
 
Vec3D operator+ (const Vec3D &a) const
 Adds another vector. More...
 
Vec3D operator- (const Vec3D a) const
 Binary vector subtraction. More...
 
bool operator== (const Vec3D &a) const
 
Vec3D multiplyElementwise (const Vec3D &a) const
 
Vec3D divideElementwise (const Vec3D &a) const
 
Vec3D signedSquare () const
 
Vec3D operator* (const Mdouble a) const
 Multiplies by a scalar. More...
 
Vec3D operator/ (Mdouble a) const
 Divides by a scalar. More...
 
Vec3Doperator+= (const Vec3D &a)
 Adds another vector. More...
 
bool operator>= (const Vec3D &a) const
 Checks if all coordinates satisfy this>=a. More...
 
bool operator< (const Vec3D &a) const
 
Vec3Doperator-= (const Vec3D &a)
 Subtracts another vector. More...
 
Vec3Doperator*= (Mdouble a)
 Multiplies by a scalar. More...
 
Vec3Doperator/= (const Mdouble a)
 Divides by a scalar. More...
 
void normalise ()
 Makes this Vec3D unit length. More...
 
void setLength (Mdouble length)
 Make this Vec3D a certain length. More...
 
Mdouble getLength () const
 Calculates the length of this Vec3D: \( \sqrt{a\cdot a} \). More...
 
Mdouble getLengthSquared () const
 Calculates the squared length of this Vec3D: \( a\cdot a \). More...
 
Mdouble getComponent (int index) const
 Returns the requested component of this Vec3D. More...
 
void setComponent (int index, double val)
 Sets the requested component of this Vec3D to the requested value. More...
 
Mdoublex ()
 RW reference to X. More...
 
Mdouble x () const
 RO reference to X. More...
 
Mdoubley ()
 RW reference to Y. More...
 
Mdouble y () const
 RO reference to Y. More...
 
Mdoublez ()
 RW reference to Z. More...
 
Mdouble z () const
 RO reference to Z. More...
 
void setX (Mdouble x)
 
void setY (Mdouble y)
 
void setZ (Mdouble z)
 
Mdouble getX () const
 
Mdouble getY () const
 
Mdouble getZ () const
 
void set (Mdouble x, Mdouble y, Mdouble z)
 
Mdouble getRadialCoordinateSquared () const
 Returns the square of the radial cylindrical coordinate, r^2=x^2+y^2. More...
 
Mdouble getRadialCoordinate () const
 Returns the square of the radial cylindrical coordinate, r=sqrt(x^2+y^2). More...
 
Vec3D getCylindricalCoordinates () const
 Returns the representation of this Vec3D in cylindrical coordinates. More...
 
Vec3D getFromCylindricalCoordinates () const
 Returns the representation of this Vec3D in cylindrical coordinates. More...
 
Vec3D getCylindricalTensorField (const Vec3D &position) const
 Returns this vector field at point p to cylindrical coordinates. More...
 
bool isEqualTo (const Vec3D &other, double tol) const
 Checks if the length this Vec3D is equal the length of other with a certain tolerance. More...
 

Static Public Member Functions

static Mdouble dot (const Vec3D &a, const Vec3D &b)
 Calculates the dot product of two Vec3D: \( a \cdot b\). More...
 
static Vec3D max (const Vec3D &a, const Vec3D &b)
 Calculates the pointwise maximum of two Vec3D. More...
 
static Vec3D min (const Vec3D &a, const Vec3D &b)
 Calculates the pointwise minimum of two Vec3D. More...
 
static double max (const Vec3D &a)
 Calculates the maximum coordinate of vector a. More...
 
static double min (const Vec3D &a)
 Calculates the minimum coordinate of vector a. More...
 
static Vec3D square (const Vec3D &a)
 Calculates the pointwise square of a Vec3D. More...
 
static Vec3D sqrt (const Vec3D &a)
 Calculates the pointwise square root of a Vec3D. More...
 
static Vec3D cross (const Vec3D &a, const Vec3D &b)
 Calculates the cross product of two Vec3D: \( a \times b\). More...
 
static Mdouble getDistance (const Vec3D &a, const Vec3D &b)
 Calculates the distance between two Vec3D: \( \sqrt{\left(a-b\right) \cdot \left(a-b\right)} \). More...
 
static Mdouble getDistanceSquared (const Vec3D &a, const Vec3D &b)
 Calculates the squared distance between two Vec3D: \( \left(a-b\right) \cdot \left(a-b\right) \). More...
 
static Mdouble getLength (const Vec3D &a)
 Calculates the length of a Vec3D: \( \sqrt{a\cdot a} \). More...
 
static Mdouble getLengthSquared (const Vec3D &a)
 Calculates the squared length of a Vec3D: \( a\cdot a \). More...
 
static Vec3D getUnitVector (const Vec3D &a)
 Returns a unit Vec3D based on a. More...
 

Public Attributes

Mdouble X
 the vector components More...
 
Mdouble Y
 
Mdouble Z
 

Friends

std::ostream & operator<< (std::ostream &os, const Vec3D &a)
 Adds elements to an output stream. More...
 
std::istream & operator>> (std::istream &is, Vec3D &a)
 Adds elements to an input stream. More...
 
Vec3D operator- (const Vec3D &a)
 Reverts the direction of a vector. More...
 
Vec3D operator* (Mdouble a, const Vec3D &b)
 Multiplies all elements by a scalar. More...
 

Constructor & Destructor Documentation

◆ Vec3D() [1/4]

Vec3D::Vec3D ( )
inline

◆ Vec3D() [2/4]

Vec3D::Vec3D ( const SmallVector< 3 > &  vector)

Alternative constructor, that constructs a Vec3D from a SmallVector size 3

Parameters
[in]vectorSmall vector that should be copied
34 {
35  X = vector[0];
36  Y = vector[1];
37  Z = vector[2];
38 }
Mdouble Y
Definition: Vector.h:66
Mdouble Z
Definition: Vector.h:66
Mdouble X
the vector components
Definition: Vector.h:66

References X, Y, and Z.

◆ Vec3D() [3/4]

Vec3D::Vec3D ( const Mdouble  x,
const Mdouble  y,
const Mdouble  z 
)
inline

Alternative constructor, taking the three elements as arguments.

Alternative constructor, lets you define all three elements.

Parameters
[in]xthe x-component
[in]ythe y-component
[in]zthe z-component
84  {
85  X = x;
86  Y = y;
87  Z = z;
88  }
Mdouble & y()
RW reference to Y.
Definition: Vector.h:372
Mdouble & z()
RW reference to Z.
Definition: Vector.h:384
Mdouble & x()
RW reference to X.
Definition: Vector.h:360

References X, x(), Y, y(), Z, and z().

◆ Vec3D() [4/4]

Vec3D::Vec3D ( std::array< double, 3 >  a)
inline

Defines mapping between std::array<double,3> and Vec3D

94  {
95  X = a[0];
96  Y = a[1];
97  Z = a[2];
98  }

References X, Y, and Z.

Member Function Documentation

◆ cross()

Vec3D Vec3D::cross ( const Vec3D a,
const Vec3D b 
)
static

Calculates the cross product of two Vec3D: \( a \times b\).

Calculates the cross product of two vectors NB: this is a STATIC function!

Parameters
[in]athe first vector
[in]bthe second vector
Returns
the cross product of the arguments
164 {
165  return Vec3D(a.Y * b.Z - a.Z * b.Y, a.Z * b.X - a.X * b.Z, a.X * b.Y - a.Y * b.X);
166 }
Vec3D()
constructor
Definition: Vector.h:71

References Vec3D(), X, Y, and Z.

Referenced by MarbleRun::actionsAfterTimeStep(), IntersectionOfWalls::add3PointObject(), ChuteWithHopper::addHopper(), IntersectionOfWalls::addTetra(), IntersectionOfWalls::addTetraSTL(), MeshTriangle::checkInteractions(), DPMBase::computeAllForces(), Mercury3Dclump::computeAllForces(), DPMBase::computeForcesDueToWalls(), Mercury3Dclump::computeForcesDueToWalls(), FrictionInteraction::computeFrictionForce(), MindlinRollingTorsionInteraction::computeFrictionForce(), DPMBase::computeInternalForce(), Mercury3Dclump::computeInternalForce(), ChuteWithPeriodicInflow::computeInternalForces(), IntersectionOfWalls::createOpenPrism(), IntersectionOfWalls::createPrism(), BaseParticle::getAngularMomentum(), MeshTriangle::getBaricentricWeight(), HorizontalBaseScrew::getDistanceAndNormal(), BasicIntersectionOfWalls::getDistanceAndNormal(), BaseWall::getInteractionWith(), TriangleMeshWall::getInteractionWith(), ParticleParticleCollision::getRelativeVelocity(), WallParticleCollision::getRelativeVelocity(), Membrane::Edge::getSineHalfTheta(), BaseInteractable::getVelocityAtContact(), TriangleMeshWall::getVolumeTetrahedron(), InfiniteWall::InfiniteWall(), TriangleWall::isInsideTriangle(), main(), WearableTriangulatedWall::processDebris(), ClumpParticle::rotatePrincipalDirections(), AngledPeriodicBoundary::set(), TriangulatedWall::setNormalsAndNeighbours(), IntersectionOfWalls::setPointsAndLines(), AxisymmetricHopper::setupInitialConditions(), WearableTriangleMeshWall::storeDebris(), WearableTriangulatedWall::storeDebris(), UniformRandomPDs(), ClumpParticle::updatePebblesVelPos(), MeshTriangle::updateVertexAndNormal(), TriangleWall::updateVertexAndNormal(), and WallVTKWriter::writeVTKSurfaceArea().

◆ divideElementwise()

Vec3D Vec3D::divideElementwise ( const Vec3D a) const
inline
152  {
153  return Vec3D(X/a.X, Y/a.Y, Z/a.Z);
154  }

References Vec3D(), X, Y, and Z.

◆ dot()

Mdouble Vec3D::dot ( const Vec3D a,
const Vec3D b 
)
static

Calculates the dot product of two Vec3D: \( a \cdot b\).

Calculates the dot product of two vectors. NB: this is a STATIC function!

Parameters
[in]athe first vector
[in]bthe second vector
Returns
the resulting scalar
77 {
78  return a.X * b.X + a.Y * b.Y + a.Z * b.Z;
79 }

References X, Y, and Z.

Referenced by IntersectionOfWalls::addTetraSTL(), Membrane::Edge::applyBendForce(), Membrane::Edge::applyStretchForce(), Membrane::Edge::calculateUPre(), ConstantMassFlowMaserBoundary::checkBoundaryAfterParticleMoved(), CGHandler::computeContactPoints(), FrictionInteraction::computeFrictionForce(), MindlinInteraction::computeFrictionForce(), MindlinRollingTorsionInteraction::computeFrictionForce(), ChuteWithPeriodicInflow::computeInternalForces(), HertzianViscoelasticInteraction::computeNormalForce(), LinearPlasticViscoelasticInteraction::computeNormalForce(), LinearViscoelasticInteraction::computeNormalForce(), SinterInteraction::computeNormalForce(), SinterLinInteraction::computeNormalForce(), SPHInteraction::computeNormalForce(), HertzianSinterInteraction::computeSinterForce(), SlidingFrictionInteraction::computeSlidingSpring(), SlidingFrictionInteraction::computeSlidingSpringSuperQuadric(), HorizontalBaseScrew::convertLimits(), AngledPeriodicBoundary::distance(), BaseInteraction::gatherContactStatistics(), MeshTriangle::getBaricentricWeight(), CGCoordinates::R::getCNormal(), ConstantMassFlowMaserBoundary::getDistance(), SubcriticalMaserBoundary::getDistance(), TriangulatedWall::Face::getDistance(), DeletionBoundary::getDistance(), FluxBoundary::getDistance(), PeriodicBoundary::getDistance(), TimeDependentPeriodicBoundary::getDistance(), NurbsSurface::getDistance(), ArcWall::getDistanceAndNormal(), AxisymmetricIntersectionOfWalls::getDistanceAndNormal(), HorizontalBaseScrew::getDistanceAndNormal(), BasicIntersectionOfWalls::getDistanceAndNormal(), Combtooth::getDistanceAndNormal(), TriangleWall::getDistanceAndNormal(), TriangulatedWall::Face::getDistanceAndNormal(), IntersectionOfWalls::getDistanceAndNormal(), SubcriticalMaserBoundaryTEST::getDistanceFromRight(), InfiniteWall::getDistanceNormalOverlapSuperquadric(), MeshTriangle::getDistanceNormalOverlapType(), DPMBase::getGravitationalEnergy(), CGCoordinates::R::getINormal(), CGCoordinates::XYZ::getINormal(), InteractionHandler::getInteraction(), BaseWall::getInteractionWith(), BaseWall::getLinePlaneIntersect(), AngledPeriodicBoundary::getOpeningAngle(), CGCoordinates::R::getPNormal(), BaseParticle::getRotationalEnergy(), ClumpParticle::getRotationalEnergy(), Membrane::Edge::getSineHalfTheta(), MindlinInteraction::getTangentialOverlap(), Membrane::getVolume(), TriangleMeshWall::getVolumeTetrahedron(), InfiniteWallWithHole::getWallDistance(), FrictionInteraction::integrate(), MindlinRollingTorsionInteraction::integrate(), PeriodicBoundary::isClosestToLeftBoundary(), TimeDependentPeriodicBoundary::isClosestToLeftBoundary(), ConstantMassFlowMaserBoundary::isClosestToRightBoundary(), SubcriticalMaserBoundary::isClosestToRightBoundary(), BaseWall::isInsideWallVTK(), main(), InfiniteWallWithHole::move_time(), SuperQuadricParticle::overlapFromContactPoint(), BaseWall::projectOntoWallVTK(), ClumpParticle::rotatePrincipalDirections(), ClumpParticle::rotateTensorOfInertia(), DeletionBoundary::set(), FluxBoundary::set(), ArcWall::set(), ConstantMassFlowMaserBoundary::set(), PeriodicBoundary::set(), SubcriticalMaserBoundary::set(), InfiniteWallWithHole::set(), TimeDependentPeriodicBoundary::set(), AngledPeriodicBoundary::set(), ChutePeriodic::set_chute_parameters(), Vreman::set_symmetric_contraction(), ContractionWithPeriodicInflow::set_symmetric_contraction(), ChuteWithPeriodicInflowAndContraction::set_symmetric_contraction(), ChuteWithContraction::set_symmetric_contraction(), ConstantMassFlowMaserBoundary::setPlanewiseShift(), PeriodicBoundary::setPlanewiseShift(), SubcriticalMaserBoundary::setPlanewiseShift(), IntersectionOfWalls::setPointsAndLines(), AngledPeriodicBoundary::shiftPosition(), AngledPeriodicBoundary::shiftPositions(), WearableTriangulatedWall::storeDebris(), ClumpParticle::updateExtraQuantities(), MindlinInteraction::updateK_t(), DPMBase::writeEneTimeStep(), LawinenBox::writeEneTimeStep(), BaseCoupling< M, O >::writeEneTimeStep(), SilbertHstop::writeToEne(), and BaseInteraction::writeToFStat().

◆ getComponent()

Mdouble Vec3D::getComponent ( int  index) const

Returns the requested component of this Vec3D.

returns the vector element belonging to the given index.

Parameters
[in]indexthe index of interest (should be 0, 1 or 2)
Returns
the value of the vector element belonging to the given index
195 {
196  switch (index)
197  {
198  case 0:
199  return X;
200  case 1:
201  return Y;
202  case 2:
203  return Z;
204  default:
205  logger(ERROR, "[Vector::getComponent] Index = %, which is too high for a 3D vector (should be 0-2).",
206  index);
207  return 0;
208  }
209 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ ERROR

References ERROR, logger, X, Y, and Z.

Referenced by MeshTriangle::checkInteractions(), Dipole::computeMultipoleExpansion(), Domain::containsParticle(), Panel::createPanels(), Domain::findNearbyBoundaries(), Panel::Panel(), Domain::setBounds(), and Domain::setRange().

◆ getCylindricalCoordinates()

Vec3D Vec3D::getCylindricalCoordinates ( ) const

Returns the representation of this Vec3D in cylindrical coordinates.

Transforms the (Cartesian) vector to cylindrical coordinates

Returns
Transformed vector
252 {
253  return Vec3D(std::sqrt(X * X + Y * Y), std::atan2(Y, X), Z);
254 }

References Vec3D(), X, Y, and Z.

◆ getCylindricalTensorField()

Vec3D Vec3D::getCylindricalTensorField ( const Vec3D p) const

Returns this vector field at point p to cylindrical coordinates.

Transforms the (Cartesian) vector to cylindrical coordinates. See https://en.wikipedia.org/wiki/Vector_fields_in_cylindrical_and_spherical_coordinates

Returns
Transformed vector
262 {
263  //define sin(A)=y/r, cos(A)=x/r
264  Mdouble r = std::sqrt(p.X * p.X + p.Y * p.Y);
265  Mdouble s = p.Y / r;
266  Mdouble c = p.X / r;
267  if (r == 0)
268  {
269  s = 0;
270  c = 1;
271  }
272  return Vec3D(X * c + Y * s, -X * s + Y * c, Z);
273 }
double Mdouble
Definition: GeneralDefine.h:34

References Vec3D(), X, Y, and Z.

Referenced by main(), CGFields::StandardFields::setCylindricalFields(), and CGFields::GradVelocityField::setCylindricalFields().

◆ getDistance()

Mdouble Vec3D::getDistance ( const Vec3D a,
const Vec3D b 
)
static

Calculates the distance between two Vec3D: \( \sqrt{\left(a-b\right) \cdot \left(a-b\right)} \).

Calculates the square of the distance (i.e. the length of the difference) between two vectors. NB: this is a STATIC function!

Parameters
[in]athe first vector
[in]bthe second vector
Returns
the square of the distance between the two arguments.

Calculates the distance (i.e. the length of the difference) between two vectors NB: this is a STATIC function!

Parameters
[in]athe first vector
[in]bthe second vector
Returns
the distance between the two arguments.
176 {
177  return std::sqrt(getDistanceSquared(a, b));
178 }
static Mdouble getDistanceSquared(const Vec3D &a, const Vec3D &b)
Calculates the squared distance between two Vec3D: .
Definition: Vector.h:311

References getDistanceSquared().

Referenced by main(), and Panel::setPanelInteractions().

◆ getDistanceSquared()

static Mdouble Vec3D::getDistanceSquared ( const Vec3D a,
const Vec3D b 
)
inlinestatic

Calculates the squared distance between two Vec3D: \( \left(a-b\right) \cdot \left(a-b\right) \).

311  {
312  const double X = a.X-b.X;
313  const double Y = a.Y-b.Y;
314  const double Z = a.Z-b.Z;
315  return (X * X + Y * Y + Z * Z);
316  //return getLengthSquared(a - b);
317  }

References X, Y, and Z.

Referenced by Mercury3Dclump::checkClumpForInteraction(), Mercury3Dclump::checkClumpForInteractionPeriodic(), DPMBase::checkParticleForInteractionLocal(), ChuteWithPeriodicInflow::computeInternalForces(), getDistance(), Mercury3D::hGridFindContactsWithTargetCell(), and BaseParticle::isInContactWith().

◆ getFromCylindricalCoordinates()

Vec3D Vec3D::getFromCylindricalCoordinates ( ) const

Returns the representation of this Vec3D in cylindrical coordinates.

Transforms the (cylindrical) vector to cartesian coordinates

Returns
Transformed vector
Todo:
280 {
282  return Vec3D(X * std::cos(Y), X * std::sin(Y), Z);
283 }
Mdouble cos(Mdouble x)
Definition: ExtendedMath.cc:64
Mdouble sin(Mdouble x)
Definition: ExtendedMath.cc:44

References mathsFunc::cos(), mathsFunc::sin(), Vec3D(), X, Y, and Z.

◆ getLength() [1/2]

Mdouble Vec3D::getLength ( ) const

Calculates the length of this Vec3D: \( \sqrt{a\cdot a} \).

Calculates the length of this vector

Returns
the (scalar) length of this vector
321 {
322  return std::sqrt(getLengthSquared());
323 }
Mdouble getLengthSquared() const
Calculates the squared length of this Vec3D: .
Definition: Vector.cc:184

References getLengthSquared().

Referenced by CGHandler::computeContactPoints(), AngledPeriodicBoundarySecondUnitTest::computeExternalForces(), ChuteWithPeriodicInflow::computeInternalForces(), SphericalWall::getDistance(), SphericalWall::getDistanceAndNormal(), TriangulatedWall::Face::getDistanceAndNormal(), AngledPeriodicBoundary::set(), AxisymmetricHopper::setupInitialConditions(), and testCGHandler().

◆ getLength() [2/2]

Mdouble Vec3D::getLength ( const Vec3D a)
static

Calculates the length of a Vec3D: \( \sqrt{a\cdot a} \).

Calculates the length of a given vector NB: this is a STATIC function!

Parameters
[in]avector to be measured.
Returns
length of the argument.
332 {
333  return a.getLength();
334 }
static Mdouble getLength(const Vec3D &a)
Calculates the length of a Vec3D: .
Definition: Vector.cc:331

References getLength().

Referenced by DrivenParticleClass::actionsAfterTimeStep(), NautaMixer::addScrew(), Membrane::Edge::applyBendForce(), BaseCluster::applyCentralForce(), Membrane::Edge::applyStretchForce(), BondedInteraction::bondInPlace(), Membrane::buildMesh(), Combtooth::Combtooth(), Mercury3Dclump::computeAllForces(), CGHandler::computeContactPoints(), MindlinInteraction::computeFrictionForce(), MindlinRollingTorsionInteraction::computeFrictionForce(), ChuteWithPeriodicInflow::computeInternalForces(), ScrewsymmetricIntersectionOfWalls::computeNormalRadialDeltaN(), SlidingFrictionInteraction::computeSlidingSpringSuperQuadric(), BaseInteraction::gatherContactStatistics(), AxisymmetricIntersectionOfWalls::getDistanceAndNormal(), HorizontalBaseScrew::getDistanceAndNormal(), Coil::getDistanceAndNormal(), Combtooth::getDistanceAndNormal(), CylindricalWall::getDistanceAndNormal(), InfiniteWallWithHole::getDistanceAndNormal(), ParabolaChute::getDistanceAndNormal(), SineWall::getDistanceAndNormal(), TriangleWall::getDistanceAndNormal(), VChute::getDistanceAndNormal(), LevelSetWall::getDistanceAndNormalLabCoordinates(), MeshTriangle::getDistanceNormalOverlapType(), getLength(), CGCoordinates::RZ::getLength(), TimeDependentPeriodicBoundary::getPlanewiseShift(), ClumpParticle::getRotationalEnergy(), MindlinInteraction::getTangentialOverlap(), SlidingFrictionInteraction::getTangentialOverlap(), main(), BaseCluster::makeDataAnalysis(), FileReader::read(), ClumpParticle::rotatePrincipalDirections(), Combtooth::set(), Chute::setChuteAngle(), WearableTriangleMeshWall::storeDebris(), ClumpParticle::updateExtraQuantities(), MindlinInteraction::updateK_t(), MeshTriangle::updateVertexAndNormal(), SingleParticle< SpeciesType >::writeEneTimeStep(), BaseInteraction::writeToFStat(), and WallVTKWriter::writeVTKSurfaceArea().

◆ getLengthSquared() [1/2]

◆ getLengthSquared() [2/2]

static Mdouble Vec3D::getLengthSquared ( const Vec3D a)
inlinestatic

Calculates the squared length of a Vec3D: \( a\cdot a \).

Calculates the square of the length of a given vector. NB: this is a STATIC function!

Parameters
[in]athe vector.
Returns
the square of the length of the argument.
333  {
334  return (a.X * a.X + a.Y * a.Y + a.Z * a.Z);
335  }

References X, Y, and Z.

Referenced by VolumeCoupling::actionsAfterTimeStep(), LawinenBox::actionsBeforeTimeStep(), IntersectionOfWalls::add3PointObject(), ChuteWithHopper::addHopper(), IntersectionOfWalls::addPlate(), IntersectionOfWalls::addTetra(), IntersectionOfWalls::addTetraSTL(), PeriodicBoundary::checkBoundaryAfterParticlesMove(), FrictionInteraction::computeFrictionForce(), MindlinRollingTorsionInteraction::computeFrictionForce(), SlidingFrictionInteraction::computeFrictionForce(), ChuteWithPeriodicInflow::computeInternalForces(), NurbsSurface::getDistance(), IntersectionOfWalls::getDistanceAndNormal(), SimpleDrumSuperquadrics::getDistanceNormalOverlapSuperquadric(), FrictionInteraction::getElasticEnergy(), MindlinInteraction::getElasticEnergy(), MindlinRollingTorsionInteraction::getElasticEnergy(), SlidingFrictionInteraction::getElasticEnergy(), BaseParticle::getKineticEnergy(), ClumpParticle::getRotationalEnergy(), getUnitVector(), InfiniteWall::InfiniteWall(), BaseParticle::integrateBeforeForceComputation(), ClumpParticle::integrateBeforeForceComputation(), LawinenBox::printTime(), setDistance(), DPMBase::setMeanVelocityAndKineticEnergy(), Domain::updateParticlePosition(), and PeriodicBoundaryHandler::updateParticles().

◆ getRadialCoordinate()

Mdouble Vec3D::getRadialCoordinate ( ) const

Returns the square of the radial cylindrical coordinate, r=sqrt(x^2+y^2).

243 {
244  return std::sqrt(X * X + Y * Y);
245 }

References X, and Y.

◆ getRadialCoordinateSquared()

Mdouble Vec3D::getRadialCoordinateSquared ( ) const

Returns the square of the radial cylindrical coordinate, r^2=x^2+y^2.

238 {
239  return X * X + Y * Y;
240 }

References X, and Y.

◆ getUnitVector()

Vec3D Vec3D::getUnitVector ( const Vec3D a)
static

Returns a unit Vec3D based on a.

Calculates the unit vector of a given vector (unless it is a vector with zero length; in that case it returns a 3D vector with each element equal to zero). NB: this is a STATIC function!

Parameters
[in]athe vector of interest
Returns
unit vector in the direction of the argument (unless the argument has length zero; in that case a zero-vector).
346 {
347  Mdouble Length2 = a.getLengthSquared();
348  if (Length2 != 0.0)
349  return a / std::sqrt(Length2);
350  else
351  return Vec3D(0, 0, 0);
352 }
static Mdouble getLengthSquared(const Vec3D &a)
Calculates the squared length of a Vec3D: .
Definition: Vector.h:332

References getLengthSquared(), and Vec3D().

Referenced by IntersectionOfWalls::createOpenPrism(), IntersectionOfWalls::createPrism(), MindlinInteraction::getTangentialOverlap(), WearableTriangulatedWall::processDebris(), TriangulatedWall::setNormalsAndNeighbours(), CGHandlerSelfTest::setupInitialConditions(), and MindlinInteraction::updateK_t().

◆ getX()

◆ getY()

◆ getZ()

◆ isEqualTo()

bool Vec3D::isEqualTo ( const Vec3D other,
double  tol 
) const

Checks if the length this Vec3D is equal the length of other with a certain tolerance.

Checks if the length of the vector is equal to the one given in the first argument (other), with a tolerance given in the second argument (tol).

Parameters
[in]otherthe 3D vector to check against
[in]tolthe tolerance
Returns
returns TRUE if the difference between the lengths of this vector and that given in the first argument (other) is smaller than the given tolerance.
295 {
296  if ((Vec3D::getLengthSquared(*this - other)) <= tol * tol)
297  {
298  return true;
299  }
300  else
301  {
302  return false;
303  }
304 }

References getLengthSquared().

Referenced by TriangleMeshWall::addToMesh(), STLTriangle::isEqualTo(), and main().

◆ isNaN()

bool Vec3D::isNaN ( ) const

Checks if ALL elements are zero.

Checks if ALL elements are zero

Returns
TRUE if ALL elements are zero
65 {
66  return std::isnan(X) || std::isnan(Y) || std::isnan(Z);
67 }

References X, Y, and Z.

◆ isZero()

bool Vec3D::isZero ( ) const
inline

Checks if ALL elements are zero.

114  { return X == 0.0 && Y == 0.0 && Z == 0.0; }

References X, Y, and Z.

Referenced by MeshTriangle::rotate(), TriangleWall::rotate(), and ArcWallUnitTest::setupInitialConditions().

◆ max() [1/2]

static double Vec3D::max ( const Vec3D a)
inlinestatic

Calculates the maximum coordinate of vector a.

265 {return std::max(std::max(a.X,a.Y),a.Z);}

References X, Y, and Z.

◆ max() [2/2]

Vec3D Vec3D::max ( const Vec3D a,
const Vec3D b 
)
static

Calculates the pointwise maximum of two Vec3D.

Calculates the pointwise maximum of two vectors. NB: this is a STATIC function!

Parameters
[in]athe first vector
[in]bthe second vector
Returns
The resulting vector, in which each element is the maximum of the equivalent elements of the arguments
90 {
91  return Vec3D(std::max(a.X, b.X), std::max(a.Y, b.Y), std::max(a.Z, b.Z));
92 }

References Vec3D(), X, Y, and Z.

Referenced by statistics_while_running< T >::auto_set_domain(), FileReader::read(), TriangleMeshWall::updateBoundingBox(), TriangleMeshWall::updateBoundingBoxGlobal(), MeshTriangle::updateVertexAndNormal(), and TriangleWall::updateVertexAndNormal().

◆ min() [1/2]

static double Vec3D::min ( const Vec3D a)
inlinestatic

Calculates the minimum coordinate of vector a.

270 {return std::min(std::min(a.X,a.Y),a.Z);}

References X, Y, and Z.

◆ min() [2/2]

Vec3D Vec3D::min ( const Vec3D a,
const Vec3D b 
)
static

Calculates the pointwise minimum of two Vec3D.

Calculates the pointwise minimum of two vectors. NB: this is a STATIC function!

Parameters
[in]athe first vector
[in]bthe second vector
Returns
The resulting vector, in which each element is the minimum of the equivalent elements of the arguments
103 {
104  return Vec3D(std::min(a.X, b.X), std::min(a.Y, b.Y), std::min(a.Z, b.Z));
105 }

References Vec3D(), X, Y, and Z.

Referenced by BaseWall::addParticlesAtWall(), statistics_while_running< T >::auto_set_domain(), FileReader::read(), TriangleMeshWall::updateBoundingBox(), TriangleMeshWall::updateBoundingBoxGlobal(), MeshTriangle::updateVertexAndNormal(), and TriangleWall::updateVertexAndNormal().

◆ multiplyElementwise()

Vec3D Vec3D::multiplyElementwise ( const Vec3D a) const
inline
148  {
149  return Vec3D(X*a.X, Y*a.Y, Z*a.Z);
150  }

References Vec3D(), X, Y, and Z.

Referenced by BaseWall::setForceControl(), and BaseWall::setVelocityControl().

◆ normalise()

◆ operator*()

Vec3D Vec3D::operator* ( const Mdouble  a) const
inline

Multiplies by a scalar.

Multiplies each element with a scalar

Parameters
[in]athe scalar to be multiplied with
Returns
the resulting vector
167  {
168  return Vec3D(X * a, Y * a, Z * a);
169  }

References Vec3D(), X, Y, and Z.

◆ operator*=()

Vec3D& Vec3D::operator*= ( Mdouble  a)
inline

Multiplies by a scalar.

Multiplies each element by a scalar

Parameters
[in]ascalar to be multiplied by
Returns
(reference to) itself, i.e. resulting vector
226  {
227  X *= a;
228  Y *= a;
229  Z *= a;
230  return *this;
231  }

References X, Y, and Z.

◆ operator+()

Vec3D Vec3D::operator+ ( const Vec3D a) const
inline

Adds another vector.

Adds vector to itself

Parameters
[in]avector to be added
Returns
resulting 3D vector
128  {
129  return Vec3D(X + a.X, Y + a.Y, Z + a.Z);
130  }

References Vec3D(), X, Y, and Z.

◆ operator+=()

Vec3D& Vec3D::operator+= ( const Vec3D a)
inline

Adds another vector.

Adds a vector to itself

Parameters
[in]avector to be added
Returns
(reference to) itself, i.e. resulting vector
188  {
189  X += a.X;
190  Y += a.Y;
191  Z += a.Z;
192  return *this;
193  }

References X, Y, and Z.

◆ operator-()

Vec3D Vec3D::operator- ( const Vec3D  a) const
inline

Binary vector subtraction.

Subtracts a vector from another vector

Parameters
[in]avector to be subtracted
Returns
resulting vector
139  {
140  return Vec3D(X - a.X, Y - a.Y, Z - a.Z);
141  };

References Vec3D(), X, Y, and Z.

◆ operator-=()

Vec3D& Vec3D::operator-= ( const Vec3D a)
inline

Subtracts another vector.

Subtracts a vector from itself

Parameters
[in]avector to be subtracted
Returns
(reference to) itself, i.e. resulting vector
213  {
214  X -= a.X;
215  Y -= a.Y;
216  Z -= a.Z;
217  return *this;
218  }

References X, Y, and Z.

◆ operator/()

Vec3D Vec3D::operator/ ( Mdouble  a) const
inline

Divides by a scalar.

Divides each element by a scalar

Parameters
[in]athe scalar to be divided by
Returns
resulting vector
177  {
178  return Vec3D(X / a, Y / a, Z / a);
179  }

References Vec3D(), X, Y, and Z.

◆ operator/=()

Vec3D& Vec3D::operator/= ( const Mdouble  a)
inline

Divides by a scalar.

Divides each element by a scalar

Parameters
[in]ascalar to be divided by
Returns
(reference to) itself, i.e. resulting vector
240  {
241  X /= a;
242  Y /= a;
243  Z /= a;
244  return *this;
245  }

References X, Y, and Z.

◆ operator<()

bool Vec3D::operator< ( const Vec3D a) const
inline
202  {
203  return X<a.X && Y<a.Y && Z<a.Z;
204  }

References X, Y, and Z.

◆ operator==()

bool Vec3D::operator== ( const Vec3D a) const
inline
144  {
145  return X == a.X and Y == a.Y and Z == a.Z;
146  };

References X, Y, and Z.

◆ operator>=()

bool Vec3D::operator>= ( const Vec3D a) const
inline

Checks if all coordinates satisfy this>=a.

198  {
199  return X>=a.X && Y>=a.Y && Z>=a.Z;
200  }

References X, Y, and Z.

◆ set()

void Vec3D::set ( Mdouble  x,
Mdouble  y,
Mdouble  z 
)
inline
412  {
413  X = x;
414  Y = y;
415  Z = z;
416  }

References X, x(), Y, y(), Z, and z().

◆ setComponent()

void Vec3D::setComponent ( int  index,
double  val 
)

Sets the requested component of this Vec3D to the requested value.

Sets the element of the vector belonging to the first argument (index) to the value given in the second argument (val).

Parameters
[in]indexindex of element of interest,
[in]valvalue to be set
218 {
219  switch (index)
220  {
221  case 0:
222  X = val;
223  break;
224  case 1:
225  Y = val;
226  break;
227  case 2:
228  Z = val;
229  break;
230  default:
231  logger(ERROR, "[Vector::setComponent] Index = %, which is too high for a 3D vector (should be 0-2).",
232  index);
233  }
234 }

References ERROR, logger, X, Y, and Z.

Referenced by ScaleCoupling< M, O >::computeCouplingForce(), Panel::createPanels(), Panel::Panel(), and Domain::setBounds().

◆ setLength()

void Vec3D::setLength ( Mdouble  length)

Make this Vec3D a certain length.

Sets the length of the vector to a given scalar (while maintaining the direction).

Parameters
[in]lengththe length to be set
139 {
140  this->normalise();
141  *this *= length;
142 }
void normalise()
Makes this Vec3D unit length.
Definition: Vector.cc:123

References normalise().

◆ setNaN()

void Vec3D::setNaN ( )

Sets all elements to NaN.

Sets each element to zero.

54 {
55  X = constants::NaN;
56  Y = constants::NaN;
57  Z = constants::NaN;
58 }
const Mdouble NaN
Definition: GeneralDefine.h:43

References constants::NaN, X, Y, and Z.

Referenced by BaseInteraction::BaseInteraction().

◆ setX()

void Vec3D::setX ( Mdouble  x)
inline
394  { X = x; }

References X, and x().

Referenced by SCoupling< M, O >::createDPMWallsFromFiniteElems().

◆ setY()

void Vec3D::setY ( Mdouble  y)
inline
397  { Y = y; }

References Y, and y().

Referenced by SCoupling< M, O >::createDPMWallsFromFiniteElems().

◆ setZ()

void Vec3D::setZ ( Mdouble  z)
inline
400  { Z = z; }

References Z, and z().

Referenced by SCoupling< M, O >::createDPMWallsFromFiniteElems().

◆ setZero()

◆ signedSquare()

Vec3D Vec3D::signedSquare ( ) const
inline
156  {
157  return Vec3D(fabs(X)*X, fabs(Y)*Y, fabs(Z)*Z);
158  }

References Vec3D(), X, Y, and Z.

◆ sqrt()

Vec3D Vec3D::sqrt ( const Vec3D a)
static

Calculates the pointwise square root of a Vec3D.

Calculates the pointwise square root of a given vector. NB: this is a STATIC function!

Parameters
[in]athe vector to be pointwise square rooted
Returns
the resulting vector, of which each element is the square root of the equivalent element of the argument.
152 {
153  return Vec3D(std::sqrt(a.X), std::sqrt(a.Y), std::sqrt(a.Z));
154 }

References Vec3D(), X, Y, and Z.

◆ square()

Vec3D Vec3D::square ( const Vec3D a)
static

Calculates the pointwise square of a Vec3D.

Calculates the pointwise square of the vector. NB: this is a STATIC function!

Parameters
[in]athe vector to be squared.
Returns
the resulting vector, of which each element is the square of the equivalent element of the argument.
115 {
116  return Vec3D(a.X * a.X, a.Y * a.Y, a.Z * a.Z);
117 }

References Vec3D(), X, Y, and Z.

Referenced by CGFields::GradVelocityField::getSquared(), and CGFields::StandardFields::getSquared().

◆ x() [1/2]

◆ x() [2/2]

Mdouble Vec3D::x ( ) const
inline

RO reference to X.

367  { return X; }

References X.

◆ y() [1/2]

◆ y() [2/2]

Mdouble Vec3D::y ( ) const
inline

RO reference to Y.

379  { return Y; }

References Y.

◆ z() [1/2]

◆ z() [2/2]

Mdouble Vec3D::z ( ) const
inline

RO reference to Z.

391  { return Z; }

References Z.

Friends And Related Function Documentation

◆ operator*

Vec3D operator* ( Mdouble  a,
const Vec3D b 
)
friend

Multiplies all elements by a scalar.

Multiplies each element of a given vector (b) by a given scalar (a). NB: this is a global function and a friend of the Vec3D class. Gets called when a scalar multiplication of the form (Mdouble) * (Vec3D) is performed.

Parameters
[in]athe scalar
[in]bthe vector
Returns
the resulting vector
479  {
480  return Vec3D(b.X * a, b.Y * a, b.Z * a);
481  }

◆ operator-

Vec3D operator- ( const Vec3D a)
friend

Reverts the direction of a vector.

466  {
467  return Vec3D(-a.X, -a.Y, -a.Z);
468  }

◆ operator<<

std::ostream& operator<< ( std::ostream &  os,
const Vec3D a 
)
friend

Adds elements to an output stream.

Adds all elements of the vector to an output stream. NB: this is a global function and a friend of the Vec3D class!

Parameters
[in]osthe output stream,
[in]aThe vector of interest
Returns
the output stream with vector elements added
362 {
363  os << a.X << ' ' << a.Y << ' ' << a.Z;
364  return os;
365 }

◆ operator>>

std::istream& operator>> ( std::istream &  is,
Vec3D a 
)
friend

Adds elements to an input stream.

Reads all elements of a given vector from an input stream. NB: this is a global function and a friend of the Vec3D class!

Parameters
[in,out]isthe input stream
[in,out]athe vector to be read in
Returns
the input stream from which the vector elements were read
375 {
376  //TW: clearing the stream avoids the nasty problem that the failbit is set to true if numbers below DBL_MIN=1e-308 are read.
377  is >> a.X;
378  if (is.fail()) {
379  logger(VERBOSE,"Failed to read x-value %", a.X);
380  //a.X = 0;
381  is.clear();
382  }
383  is >> a.Y;
384  if (is.fail()) {
385  logger(VERBOSE,"Failed to read y-value %", a.Y);
386  //a.Y = 0;
387  is.clear();
388  }
389  is >> a.Z;
390  if (is.fail()) {
391  logger(VERBOSE,"Failed to read z-value %", a.Z);
392  //a.Z = 0;
393  is.clear();
394  }
395  return is;
396 }
@ VERBOSE

Member Data Documentation

◆ X

Mdouble Vec3D::X

the vector components

Referenced by ScaleCoupledBeam::actionsAfterSolve(), HertzianBSHPInteractionTwoParticleElasticCollision::actionsAfterSolve(), SphericalSuperQuadricCollision::actionsAfterSolve(), SpeciesTest::actionsAfterSolve(), T_protectiveWall::actionsAfterTimeStep(), ContactDetectionNormalSpheresTest::actionsAfterTimeStep(), ContactDetectionRotatedSpheresTest::actionsAfterTimeStep(), protectiveWall::actionsAfterTimeStep(), DrivenParticleClass::actionsAfterTimeStep(), LawinenBox::actionsBeforeTimeStep(), SmoothChute::actionsBeforeTimeStep(), AngleOfRepose::actionsBeforeTimeStep(), Chutebelt::actionsBeforeTimeStep(), StressStrainControlBoundary::activateStrainRateControl(), ChutePeriodic::add_particles(), NautaMixer::addBaseWall(), NautaMixer::addConeWall(), ChuteWithHopper::addHopper(), NautaMixer::addParticles(), NautaMixer::addParticlesAtWall(), BaseWall::addParticlesAtWall(), ClumpParticle::angularAccelerateClumpIterative(), Quaternion::angularDisplacementTimeDerivative(), Quaternion::angularVelocityBodyFixedFrameToAngularDisplacement(), statistics_while_running< T >::auto_set_domain(), DeletionBoundary::checkBoundaryAfterParticleMoved(), HeaterBoundary::checkBoundaryAfterParticleMoved(), VolumeCoupling::checkParticlesInFiniteElems(), ChuteWithContraction::ChuteWithContraction(), ChuteWithPeriodicInflow::cleanChute(), ChuteWithContraction::cleanChute(), Chute::cleanChute(), ScrewsymmetricIntersectionOfWalls::computeDeltaZ(), ScaleCoupling< M, O >::computeExternalForces(), SuperQuadricParticle::computeHessianLabFixed(), Mercury2D::computeInternalForces(), Mercury3D::computeInternalForces(), ChuteWithPeriodicInflow::computeInternalForces(), BaseCluster::computeInternalStructure(), DPM::computeLocalCGHGrid(), DPM::computeLocalVolumeFractionHGrid(), SuperQuadricParticle::computeMass(), ScrewsymmetricIntersectionOfWalls::computeNormalRadialDeltaN(), SuperQuadricParticle::computeResidualContactDetection(), SuperQuadricParticle::computeShape(), SuperQuadricParticle::computeShapeGradientLabFixed(), Mercury3D::computeWallForces(), AxisymmetricIntersectionOfWalls::convertLimits(), HorizontalBaseScrew::convertLimits(), ScrewsymmetricIntersectionOfWalls::convertLimits(), Funnel::create_funnel(), LawinenBox::create_inflow_particle(), Funnel::create_inflow_particle(), FlowRule::create_inflow_particle(), SilbertPeriodic::create_inflow_particle(), Chute::createBottom(), CircularPeriodicBoundary::createPeriodicParticle(), InfiniteWall::createVTK(), Matrix3D::cross(), cross(), divideElementwise(), dot(), Matrix3D::dyadic(), ChuteWithWedge::extendBottom(), SubcriticalMaserBoundaryTEST::extendBottom(), MembraneDemo::fixMembraneEdges(), BidisperseCubeInsertionBoundary::generateParticle(), BaseWall::getAxis(), MeshTriangle::getBaricentricWeight(), CGCoordinates::RZ::getCNormal(), CGCoordinates::X::getCNormal(), CGCoordinates::XY::getCNormal(), CGCoordinates::XZ::getCNormal(), getComponent(), SuperQuadricParticle::getContactPointPlanB(), getCylindricalCoordinates(), Matrix3D::getCylindricalTensorField(), getCylindricalTensorField(), BaseParticle::getDisplacement2(), HeaterBoundary::getDistance(), CubeDeletionBoundary::getDistance(), Quaternion::getDistance(), AxisymmetricIntersectionOfWalls::getDistanceAndNormal(), HorizontalBaseScrew::getDistanceAndNormal(), HorizontalScrew::getDistanceAndNormal(), ScrewsymmetricIntersectionOfWalls::getDistanceAndNormal(), BasicIntersectionOfWalls::getDistanceAndNormal(), Coil::getDistanceAndNormal(), CylindricalWall::getDistanceAndNormal(), InfiniteWallWithHole::getDistanceAndNormal(), SineWall::getDistanceAndNormal(), VChute::getDistanceAndNormal(), LevelSetWall::getDistanceAndNormalLabCoordinates(), Screw::getDistanceAndNormalLabCoordinates(), getDistanceSquared(), CGCoordinates::RZ::getDistanceSquared(), CGCoordinates::X::getDistanceSquared(), CGCoordinates::XY::getDistanceSquared(), CGCoordinates::XZ::getDistanceSquared(), CGCoordinates::BaseCoordinates::getDomainVolume(), CGCoordinates::R::getDomainVolume(), SCoupling< M, O >::getElementBoundingBox(), getFromCylindricalCoordinates(), InfiniteWall::getFurthestPointSuperQuadric(), InfiniteWallWithHole::getHoleDistance(), LeesEdwardsBoundary::getHorizontalDistance(), ShearBoxBoundary::getHorizontalDistance(), WallSpecies::getInfo(), SuperQuadricParticle::getInitialGuessForContact(), CGCoordinates::RZ::getINormal(), CGCoordinates::X::getINormal(), CGCoordinates::XY::getINormal(), CGCoordinates::XZ::getINormal(), SuperQuadricParticle::getInteractionWithSuperQuad(), SuperQuadricParticle::getJacobianOfContactDetectionObjective(), ClumpParticle::getKineticEnergy(), CGCoordinates::R::getLength(), CGCoordinates::X::getLength(), CGCoordinates::XY::getLength(), CGCoordinates::XYZ::getLength(), CGCoordinates::XZ::getLength(), getLengthSquared(), ParticleHandler::getMassTimesPosition(), getMPISum(), DomainHandler::getParticleDomainGlobalIndex(), CGCoordinates::RZ::getPNormal(), CGCoordinates::X::getPNormal(), CGCoordinates::XY::getPNormal(), CGCoordinates::XZ::getPNormal(), getRadialCoordinate(), getRadialCoordinateSquared(), ParticleParticleCollision::getRelativeVelocityComponents(), WallParticleCollision::getRelativeVelocityComponents(), CGCoordinates::RZ::getTangentialSquared(), CGCoordinates::XY::getTangentialSquared(), CGCoordinates::XZ::getTangentialSquared(), HeaterBoundary::getVolume(), SuperQuadricParticle::getVolume(), CGCoordinates::O::getVolumeOfAveragedDimensions(), CGCoordinates::Y::getVolumeOfAveragedDimensions(), CGCoordinates::YZ::getVolumeOfAveragedDimensions(), CGCoordinates::Z::getVolumeOfAveragedDimensions(), WearableNurbsWall::getVolumeUnderSurface(), WearableNurbsWall::getVolumeUnderSurfaceX(), getX(), DPMBase::getXCenter(), Mercury2D::hGridFindParticleContacts(), Mercury3D::hGridFindParticleContacts(), Mercury2D::hGridGetInteractingParticleList(), Mercury3D::hGridGetInteractingParticleList(), Mercury2D::hGridHasParticleContacts(), Mercury3D::hGridHasParticleContacts(), Mercury2D::hGridUpdateParticle(), Mercury3D::hGridUpdateParticle(), MarbleRun::includeInDomain(), InitialConditions< SpeciesType >::InitialConditions(), HorizontalMixer::introduceParticlesAtWall(), ContactDetectionIntersectionOfWallsTest::introduceParticlesAtWall(), HorizontalMixer::introduceParticlesInDomain(), mathsFunc::isEqual(), SuperQuadricParticle::isInContactWith(), MeshTriangle::isInsideTriangle(), isNaN(), isZero(), LawinenBox::LawinenBox(), Matrix3D::ldivide(), load(), main(), max(), min(), multiplyElementwise(), operator*(), Matrix3D::operator*(), operator*=(), operator+(), operator+=(), operator-(), operator-=(), operator/(), operator/=(), operator<(), operator==(), operator>=(), LiquidMigrationMPI2Test::outputXBallsData(), LiquidMigrationSelfTest::outputXBallsData(), ChuteWithPeriodicInflow::outputXBallsDataParticlee(), SuperQuadricParticle::overlapFromContactPoint(), BaseCluster::particleInsertionSuccessful(), FixedClusterInsertionBoundary::placeParticle(), ChuteInsertionBoundary::placeParticle(), CubeInsertionBoundary::placeParticle(), PolydisperseInsertionBoundary::placeParticle(), RandomClusterInsertionBoundary::placeParticle(), ChuteWithPeriodicInflow::printTime(), RotatingDrumWet::printTime(), ForceLawsMPI2Test::printTime(), BaseCluster::printTime(), DPMBase::readNextDataFile(), WallHandler::readTriangleWall(), TriangulatedWall::readVTK(), save(), MatrixSymmetric3D::selfDyadic(), serialize(), set(), AngledPeriodicBoundary::set(), SuperQuadricParticle::setBoundingRadius(), setComponent(), DPMBase::setDomain(), Quaternion::setEuler(), CGFields::OrientationField::setFields(), BaseCG::setHX(), SuperQuadricParticle::setInertia(), setNaN(), Quaternion::setOrientationViaNormal(), ClumpParticle::setPrincipalDirections_e1(), ClumpParticle::setPrincipalDirections_e2(), ClumpParticle::setPrincipalDirections_e3(), Domain::setRange(), HeaterBoundary::setStrength2D(), HeaterBoundary::setStrength3D(), ClosedCSCWalls::setupInitialConditions(), CSCInit::setupInitialConditions(), CSCWalls::setupInitialConditions(), MercuryLogo::setupInitialConditions(), SmoothChute::setupInitialConditions(), HourGlass::setupInitialConditions(), MinimalExampleDrum::setupInitialConditions(), TimeDependentPeriodicBoundary3DSelfTest::setupInitialConditions(), AngleOfRepose::setupInitialConditions(), GetDistanceAndNormalForIntersectionOfWalls::setupInitialConditions(), GetDistanceAndNormalForScrew::setupInitialConditions(), GetDistanceAndNormalForTriangleWall::setupInitialConditions(), LeesEdwardsSelfTest::setupInitialConditions(), ParticleCreation::setupInitialConditions(), my_problem_HGRID::setupInitialConditions(), GetDistanceAndNormalForTriangleWalls::setupInitialConditions(), TriangulatedScrewSelfTest::setupInitialConditions(), TriangulatedWallSelfTest::setupInitialConditions(), DrumRot::setupInitialConditions(), RotatingDrum::setupInitialConditions(), EllipticalSuperQuadricCollision::setupInitialConditions(), Tutorial11::setupInitialConditions(), ArcWallUnitTest::setupInitialConditions(), MD_demo::setupInitialConditions(), ChuteBottom::setupInitialConditions(), BaseCG::setX(), setX(), setZero(), signedSquare(), CGCoordinates::spaceEvenly(), sqrt(), square(), FlowFrontChute::stretch(), MatrixSymmetric3D::symmetrisedDyadic(), ContactDetectionNormalSpheresTest::test(), Packing::test(), HertzContactRestitutionUnitTest::test(), UniformRandomPDs(), TriangleMeshWall::updateBoundingBoxGlobal(), ScaleCoupling< M, O >::updateCoupledElements(), StressStrainControlBoundary::updateDomainSize(), PeriodicBoundaryHandler::updateMaserParticle(), Vec3D(), SuperQuadricParticle::writeDebugMessageStep2(), SuperQuadricParticle::writeDebugMessageStep3(), DPMBase::writeEneTimeStep(), LawinenBox::writeEneTimeStep(), Drum::writeEneTimeStep(), Slide::writeEneTimeStep(), BaseCoupling< M, O >::writeEneTimeStep(), CFile::writeP4C(), GranuDrum::writeResults(), BaseCluster::writeToCdatFile(), HorizontalScrew::writeVTK(), Screw::writeVTK(), and x().

◆ Y

Mdouble Vec3D::Y

Referenced by ContactDetectionNormalSpheresTest::actionsAfterTimeStep(), ContactDetectionRotatedSpheresTest::actionsAfterTimeStep(), DrivenParticleClass::actionsAfterTimeStep(), SmoothChute::actionsBeforeTimeStep(), AngleOfRepose::actionsBeforeTimeStep(), Slide::actionsBeforeTimeStep(), StressStrainControlBoundary::activateStrainRateControl(), ChuteWithHopper::addHopper(), NautaMixer::addParticles(), NautaMixer::addParticlesAtWall(), BaseWall::addParticlesAtWall(), ClumpParticle::angularAccelerateClumpIterative(), Quaternion::angularDisplacementTimeDerivative(), Quaternion::angularVelocityBodyFixedFrameToAngularDisplacement(), statistics_while_running< T >::auto_set_domain(), DeletionBoundary::checkBoundaryAfterParticleMoved(), HeaterBoundary::checkBoundaryAfterParticleMoved(), VolumeCoupling::checkParticlesInFiniteElems(), ScrewsymmetricIntersectionOfWalls::computeDeltaZ(), ScaleCoupling< M, O >::computeExternalForces(), SuperQuadricParticle::computeHessianLabFixed(), Mercury2D::computeInternalForces(), Mercury3D::computeInternalForces(), BaseCluster::computeInternalStructure(), DPM::computeLocalCGHGrid(), DPM::computeLocalVolumeFractionHGrid(), SuperQuadricParticle::computeMass(), ScrewsymmetricIntersectionOfWalls::computeNormalRadialDeltaN(), SuperQuadricParticle::computeResidualContactDetection(), SuperQuadricParticle::computeShape(), SuperQuadricParticle::computeShapeGradientLabFixed(), Mercury3D::computeWallForces(), AxisymmetricIntersectionOfWalls::convertLimits(), HorizontalBaseScrew::convertLimits(), ScrewsymmetricIntersectionOfWalls::convertLimits(), Funnel::create_funnel(), LawinenBox::create_inflow_particle(), SegregationWithHopper::create_inflow_particle(), Funnel::create_inflow_particle(), FlowRule::create_inflow_particle(), SilbertPeriodic::create_inflow_particle(), Chute::createBottom(), CircularPeriodicBoundary::createPeriodicParticle(), InfiniteWall::createVTK(), Matrix3D::cross(), cross(), divideElementwise(), dot(), Matrix3D::dyadic(), SubcriticalMaserBoundaryTEST::extendBottom(), MembraneDemo::fixMembraneEdges(), MembraneSelfTest::fixMembraneEdges(), BidisperseCubeInsertionBoundary::generateParticle(), BaseWall::getAxis(), MeshTriangle::getBaricentricWeight(), CGCoordinates::XY::getCNormal(), CGCoordinates::Y::getCNormal(), CGCoordinates::YZ::getCNormal(), getComponent(), SuperQuadricParticle::getContactPointPlanB(), getCylindricalCoordinates(), Matrix3D::getCylindricalTensorField(), getCylindricalTensorField(), BaseParticle::getDisplacement2(), HeaterBoundary::getDistance(), CubeDeletionBoundary::getDistance(), Quaternion::getDistance(), HorizontalScrew::getDistanceAndNormal(), ScrewsymmetricIntersectionOfWalls::getDistanceAndNormal(), BasicIntersectionOfWalls::getDistanceAndNormal(), Coil::getDistanceAndNormal(), CylindricalWall::getDistanceAndNormal(), InfiniteWallWithHole::getDistanceAndNormal(), ParabolaChute::getDistanceAndNormal(), SimpleDrumSuperquadrics::getDistanceAndNormal(), SineWall::getDistanceAndNormal(), VChute::getDistanceAndNormal(), LevelSetWall::getDistanceAndNormalLabCoordinates(), Screw::getDistanceAndNormalLabCoordinates(), SimpleDrumSuperquadrics::getDistanceNormalOverlapSuperquadric(), getDistanceSquared(), CGCoordinates::RZ::getDistanceSquared(), CGCoordinates::XY::getDistanceSquared(), CGCoordinates::Y::getDistanceSquared(), CGCoordinates::YZ::getDistanceSquared(), CGCoordinates::BaseCoordinates::getDomainVolume(), SCoupling< M, O >::getElementBoundingBox(), getFromCylindricalCoordinates(), InfiniteWall::getFurthestPointSuperQuadric(), InfiniteWallWithHole::getHoleDistance(), SuperQuadricParticle::getInitialGuessForContact(), CGCoordinates::XY::getINormal(), CGCoordinates::Y::getINormal(), CGCoordinates::YZ::getINormal(), SuperQuadricParticle::getInteractionWithSuperQuad(), SuperQuadricParticle::getJacobianOfContactDetectionObjective(), ClumpParticle::getKineticEnergy(), CGCoordinates::R::getLength(), CGCoordinates::XY::getLength(), CGCoordinates::XYZ::getLength(), CGCoordinates::Y::getLength(), CGCoordinates::YZ::getLength(), getLengthSquared(), ParticleHandler::getMassTimesPosition(), getMPISum(), DomainHandler::getParticleDomainGlobalIndex(), CGCoordinates::XY::getPNormal(), CGCoordinates::Y::getPNormal(), CGCoordinates::YZ::getPNormal(), getRadialCoordinate(), getRadialCoordinateSquared(), ParticleParticleCollision::getRelativeVelocityComponents(), WallParticleCollision::getRelativeVelocityComponents(), CGCoordinates::XY::getTangentialSquared(), CGCoordinates::YZ::getTangentialSquared(), LeesEdwardsBoundary::getVerticalDistance(), ShearBoxBoundary::getVerticalDistance(), HeaterBoundary::getVolume(), SuperQuadricParticle::getVolume(), CGCoordinates::O::getVolumeOfAveragedDimensions(), CGCoordinates::X::getVolumeOfAveragedDimensions(), CGCoordinates::XZ::getVolumeOfAveragedDimensions(), CGCoordinates::Z::getVolumeOfAveragedDimensions(), WearableNurbsWall::getVolumeUnderSurface(), WearableNurbsWall::getVolumeUnderSurfaceX(), getY(), DPMBase::getYCenter(), Mercury2D::hGridFindParticleContacts(), Mercury3D::hGridFindParticleContacts(), Mercury2D::hGridGetInteractingParticleList(), Mercury3D::hGridGetInteractingParticleList(), Mercury2D::hGridHasParticleContacts(), Mercury3D::hGridHasParticleContacts(), Mercury2D::hGridUpdateParticle(), Mercury3D::hGridUpdateParticle(), MarbleRun::includeInDomain(), InitialConditions< SpeciesType >::InitialConditions(), HorizontalMixer::introduceParticlesAtWall(), ContactDetectionIntersectionOfWallsTest::introduceParticlesAtWall(), HorizontalMixer::introduceParticlesInDomain(), mathsFunc::isEqual(), SuperQuadricParticle::isInContactWith(), MeshTriangle::isInsideTriangle(), isNaN(), isZero(), LawinenBox::LawinenBox(), Matrix3D::ldivide(), load(), main(), max(), min(), multiplyElementwise(), operator*(), Matrix3D::operator*(), operator*=(), operator+(), operator+=(), operator-(), operator-=(), operator/(), operator/=(), operator<(), operator==(), operator>=(), LiquidMigrationMPI2Test::outputXBallsData(), LiquidMigrationSelfTest::outputXBallsData(), ChuteWithPeriodicInflow::outputXBallsDataParticlee(), SuperQuadricParticle::overlapFromContactPoint(), BaseCluster::particleInsertionSuccessful(), FixedClusterInsertionBoundary::placeParticle(), ChuteInsertionBoundary::placeParticle(), CubeInsertionBoundary::placeParticle(), PolydisperseInsertionBoundary::placeParticle(), RandomClusterInsertionBoundary::placeParticle(), BaseCluster::printTime(), DPMBase::readNextDataFile(), WallHandler::readTriangleWall(), TriangulatedWall::readVTK(), save(), MatrixSymmetric3D::selfDyadic(), serialize(), set(), AngledPeriodicBoundary::set(), SuperQuadricParticle::setBoundingRadius(), setComponent(), DPMBase::setDomain(), Quaternion::setEuler(), CGFields::OrientationField::setFields(), BaseCG::setHY(), SuperQuadricParticle::setInertia(), setNaN(), Quaternion::setOrientationViaNormal(), ClumpParticle::setPrincipalDirections_e1(), ClumpParticle::setPrincipalDirections_e2(), ClumpParticle::setPrincipalDirections_e3(), Domain::setRange(), HeaterBoundary::setStrength2D(), HeaterBoundary::setStrength3D(), ClosedCSCWalls::setupInitialConditions(), CSCInit::setupInitialConditions(), CSCWalls::setupInitialConditions(), MercuryLogo::setupInitialConditions(), SmoothChute::setupInitialConditions(), HourGlass::setupInitialConditions(), MinimalExampleDrum::setupInitialConditions(), TimeDependentPeriodicBoundary3DSelfTest::setupInitialConditions(), AngleOfRepose::setupInitialConditions(), LeesEdwardsSelfTest::setupInitialConditions(), ParticleCreation::setupInitialConditions(), my_problem_HGRID::setupInitialConditions(), TriangulatedScrewSelfTest::setupInitialConditions(), TriangulatedWallSelfTest::setupInitialConditions(), RotatingDrum::setupInitialConditions(), EllipticalSuperQuadricCollision::setupInitialConditions(), Tutorial11::setupInitialConditions(), ArcWallUnitTest::setupInitialConditions(), MD_demo::setupInitialConditions(), ChuteBottom::setupInitialConditions(), BaseCG::setY(), setY(), setZero(), signedSquare(), CGCoordinates::spaceEvenly(), sqrt(), square(), MatrixSymmetric3D::symmetrisedDyadic(), Packing::test(), UniformRandomPDs(), TriangleMeshWall::updateBoundingBoxGlobal(), ScaleCoupling< M, O >::updateCoupledElements(), StressStrainControlBoundary::updateDomainSize(), Vec3D(), SuperQuadricParticle::writeDebugMessageStep2(), SuperQuadricParticle::writeDebugMessageStep3(), DPMBase::writeEneTimeStep(), Slide::writeEneTimeStep(), BaseCoupling< M, O >::writeEneTimeStep(), CFile::writeP4C(), BaseCluster::writeToCdatFile(), HorizontalScrew::writeVTK(), Screw::writeVTK(), and y().

◆ Z

Mdouble Vec3D::Z

Referenced by GranuHeap::actionsAfterTimeStep(), ClosedCSCRestart::actionsAfterTimeStep(), ClosedCSCRun::actionsAfterTimeStep(), ClosedCSCWalls::actionsAfterTimeStep(), MarbleRun::actionsAfterTimeStep(), BouncingSuperQuadric::actionsAfterTimeStep(), ContactDetectionNormalSpheresTest::actionsAfterTimeStep(), ContactDetectionRotatedSpheresTest::actionsAfterTimeStep(), DrivenParticleClass::actionsAfterTimeStep(), LawinenBox::actionsBeforeTimeStep(), SmoothChute::actionsBeforeTimeStep(), Slide::actionsBeforeTimeStep(), StressStrainControlBoundary::activateStrainRateControl(), ChutePeriodic::add_particles(), NautaMixer::addBaseWall(), NautaMixer::addConeWall(), ChuteWithHopper::addHopper(), NautaMixer::addParticles(), NautaMixer::addParticlesAtWall(), BaseWall::addParticlesAtWall(), ClumpParticle::angularAccelerateClumpIterative(), Quaternion::angularDisplacementTimeDerivative(), Quaternion::angularVelocityBodyFixedFrameToAngularDisplacement(), statistics_while_running< T >::auto_set_domain(), statistics_while_running< T >::auto_set_z(), DeletionBoundary::checkBoundaryAfterParticleMoved(), HeaterBoundary::checkBoundaryAfterParticleMoved(), VolumeCoupling::checkParticlesInFiniteElems(), Funnel::cleanChute(), ScrewsymmetricIntersectionOfWalls::computeDeltaZ(), ScaleCoupling< M, O >::computeExternalForces(), SuperQuadricParticle::computeHessianLabFixed(), Mercury3D::computeInternalForces(), BaseCluster::computeInternalStructure(), DPM::computeLocalCGHGrid(), DPM::computeLocalVolumeFractionHGrid(), SuperQuadricParticle::computeMass(), ScrewsymmetricIntersectionOfWalls::computeNormalRadialDeltaN(), SuperQuadricParticle::computeResidualContactDetection(), SuperQuadricParticle::computeShape(), SuperQuadricParticle::computeShapeGradientLabFixed(), Mercury3D::computeWallForces(), ClosedCSCWalls::continueSolve(), AxisymmetricIntersectionOfWalls::convertLimits(), HorizontalBaseScrew::convertLimits(), ScrewsymmetricIntersectionOfWalls::convertLimits(), Funnel::create_funnel(), LawinenBox::create_inflow_particle(), Funnel::create_inflow_particle(), FlowRule::create_inflow_particle(), SilbertPeriodic::create_inflow_particle(), InfiniteWall::createVTK(), Matrix3D::cross(), cross(), divideElementwise(), dot(), Matrix3D::dyadic(), SubcriticalMaserBoundaryTEST::extendBottom(), BidisperseCubeInsertionBoundary::generateParticle(), BaseWall::getAxis(), MeshTriangle::getBaricentricWeight(), CGCoordinates::RZ::getCNormal(), CGCoordinates::XZ::getCNormal(), CGCoordinates::YZ::getCNormal(), CGCoordinates::Z::getCNormal(), getComponent(), SuperQuadricParticle::getContactPointPlanB(), getCylindricalCoordinates(), getCylindricalTensorField(), BaseParticle::getDisplacement2(), HeaterBoundary::getDistance(), CubeDeletionBoundary::getDistance(), Quaternion::getDistance(), AxisymmetricIntersectionOfWalls::getDistanceAndNormal(), HorizontalBaseScrew::getDistanceAndNormal(), HorizontalScrew::getDistanceAndNormal(), ScrewsymmetricIntersectionOfWalls::getDistanceAndNormal(), BasicIntersectionOfWalls::getDistanceAndNormal(), Coil::getDistanceAndNormal(), CylindricalWall::getDistanceAndNormal(), InfiniteWallWithHole::getDistanceAndNormal(), ParabolaChute::getDistanceAndNormal(), SineWall::getDistanceAndNormal(), VChute::getDistanceAndNormal(), LevelSetWall::getDistanceAndNormalLabCoordinates(), Screw::getDistanceAndNormalLabCoordinates(), getDistanceSquared(), CGCoordinates::RZ::getDistanceSquared(), CGCoordinates::XZ::getDistanceSquared(), CGCoordinates::YZ::getDistanceSquared(), CGCoordinates::Z::getDistanceSquared(), CGCoordinates::BaseCoordinates::getDomainVolume(), CGCoordinates::R::getDomainVolume(), SCoupling< M, O >::getElementBoundingBox(), SphericalIndenter::getForceOnIndenter(), getFromCylindricalCoordinates(), InfiniteWall::getFurthestPointSuperQuadric(), SphericalIndenter::getIndenterHeight(), SphericalIndenter::getIndenterVelocity(), WallSpecies::getInfo(), SuperQuadricParticle::getInitialGuessForContact(), CGCoordinates::RZ::getINormal(), CGCoordinates::XZ::getINormal(), CGCoordinates::YZ::getINormal(), CGCoordinates::Z::getINormal(), SuperQuadricParticle::getInteractionWithSuperQuad(), SuperQuadricParticle::getJacobianOfContactDetectionObjective(), ClumpParticle::getKineticEnergy(), CGCoordinates::XYZ::getLength(), CGCoordinates::XZ::getLength(), CGCoordinates::YZ::getLength(), CGCoordinates::Z::getLength(), getLengthSquared(), ParticleHandler::getMassTimesPosition(), getMPISum(), DomainHandler::getParticleDomainGlobalIndex(), CGCoordinates::RZ::getPNormal(), CGCoordinates::XZ::getPNormal(), CGCoordinates::YZ::getPNormal(), CGCoordinates::Z::getPNormal(), CGCoordinates::RZ::getTangentialSquared(), CGCoordinates::XZ::getTangentialSquared(), CGCoordinates::YZ::getTangentialSquared(), HeaterBoundary::getVolume(), SuperQuadricParticle::getVolume(), CGCoordinates::O::getVolumeOfAveragedDimensions(), CGCoordinates::R::getVolumeOfAveragedDimensions(), CGCoordinates::X::getVolumeOfAveragedDimensions(), CGCoordinates::XY::getVolumeOfAveragedDimensions(), CGCoordinates::Y::getVolumeOfAveragedDimensions(), WearableNurbsWall::getVolumeUnderSurface(), WearableNurbsWall::getVolumeUnderSurfaceX(), getZ(), DPMBase::getZCenter(), GranuHeap::GranuHeap(), Mercury3D::hGridFindParticleContacts(), Mercury3D::hGridGetInteractingParticleList(), Mercury3D::hGridHasParticleContacts(), Mercury3D::hGridUpdateParticle(), MarbleRun::includeInDomain(), InitialConditions< SpeciesType >::InitialConditions(), HorizontalMixer::introduceParticlesAtWall(), ContactDetectionIntersectionOfWallsTest::introduceParticlesAtWall(), HorizontalMixer::introduceParticlesInDomain(), mathsFunc::isEqual(), SuperQuadricParticle::isInContactWith(), MeshTriangle::isInsideTriangle(), isNaN(), isZero(), Matrix3D::ldivide(), load(), main(), max(), min(), multiplyElementwise(), operator*(), Matrix3D::operator*(), operator*=(), operator+(), operator+=(), operator-(), operator-=(), operator/(), operator/=(), operator<(), operator==(), operator>=(), LiquidMigrationMPI2Test::outputXBallsData(), LiquidMigrationSelfTest::outputXBallsData(), ChuteWithPeriodicInflow::outputXBallsDataParticlee(), SuperQuadricParticle::overlapFromContactPoint(), BaseCluster::particleInsertionSuccessful(), FixedClusterInsertionBoundary::placeParticle(), ChuteInsertionBoundary::placeParticle(), CubeInsertionBoundary::placeParticle(), PolydisperseInsertionBoundary::placeParticle(), RandomClusterInsertionBoundary::placeParticle(), ClosedCSCRestart::printTime(), ClosedCSCRun::printTime(), ClosedCSCWalls::printTime(), RotatingDrumWet::printTime(), CubeDeletionBoundarySelfTest::printTime(), DeletionBoundarySelfTest::printTime(), BaseCluster::printTime(), DPMBase::readNextDataFile(), WallHandler::readTriangleWall(), TriangulatedWall::readVTK(), save(), ClosedCSCWalls::saveWalls(), MatrixSymmetric3D::selfDyadic(), serialize(), set(), AngledPeriodicBoundary::set(), SuperQuadricParticle::setBoundingRadius(), setComponent(), DPMBase::setDomain(), Quaternion::setEuler(), CGFields::OrientationField::setFields(), BaseCG::setHZ(), SuperQuadricParticle::setInertia(), setNaN(), Quaternion::setOrientationViaNormal(), ClumpParticle::setPrincipalDirections_e1(), ClumpParticle::setPrincipalDirections_e2(), ClumpParticle::setPrincipalDirections_e3(), Domain::setRange(), HeaterBoundary::setStrength2D(), HeaterBoundary::setStrength3D(), ClosedCSCWalls::setupInitialConditions(), CSCInit::setupInitialConditions(), CSCWalls::setupInitialConditions(), MercuryLogo::setupInitialConditions(), SmoothChute::setupInitialConditions(), HourGlass::setupInitialConditions(), MinimalExampleDrum::setupInitialConditions(), TimeDependentPeriodicBoundary3DSelfTest::setupInitialConditions(), GetDistanceAndNormalForIntersectionOfWalls::setupInitialConditions(), GetDistanceAndNormalForScrew::setupInitialConditions(), GetDistanceAndNormalForTriangleWall::setupInitialConditions(), NewtonsCradleSelfTest::setupInitialConditions(), ParticleCreation::setupInitialConditions(), GetDistanceAndNormalForTriangleWalls::setupInitialConditions(), TriangulatedScrewSelfTest::setupInitialConditions(), TriangulatedWallSelfTest::setupInitialConditions(), DrumRot::setupInitialConditions(), RotatingDrum::setupInitialConditions(), InitialConditions< SpeciesType >::setupInitialConditions(), EllipticalSuperQuadricCollision::setupInitialConditions(), Tutorial11::setupInitialConditions(), MD_demo::setupInitialConditions(), ChuteBottom::setupInitialConditions(), BaseCG::setZ(), setZ(), setZero(), ShearBoxBoundary::shiftHorizontalPosition(), signedSquare(), CGCoordinates::spaceEvenly(), sqrt(), square(), MatrixSymmetric3D::symmetrisedDyadic(), Packing::test(), UniformRandomPDs(), TriangleMeshWall::updateBoundingBoxGlobal(), ScaleCoupling< M, O >::updateCoupledElements(), StressStrainControlBoundary::updateDomainSize(), Vec3D(), SuperQuadricParticle::writeDebugMessageStep2(), SuperQuadricParticle::writeDebugMessageStep3(), DPMBase::writeEneTimeStep(), LawinenBox::writeEneTimeStep(), Drum::writeEneTimeStep(), Penetration::writeEneTimeStep(), SingleParticle< SpeciesType >::writeEneTimeStep(), BaseCoupling< M, O >::writeEneTimeStep(), CFile::writeP4C(), GranuDrum::writeResults(), BaseCluster::writeToCdatFile(), HorizontalBaseScrew::writeVTK(), HorizontalScrew::writeVTK(), Screw::writeVTK(), WearableNurbsWall::writeWallDetailsVTK(), and z().


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