MercuryDPM
Trunk
|
Class that implements particles which store both temperature/heat capacity and liquid content which is adapted for the CFD-DEM studies. More...
#include <HeatFluidCoupledParticle.h>
Public Member Functions | |
HeatFluidCoupledParticle () | |
HeatFluidCoupledParticle (const HeatFluidCoupledParticle &p) | |
HeatFluidCoupledParticle copy constructor, which accepts as input a reference to a HeatFluidCoupledParticle. It creates a copy of this HeatFluidCoupledParticle and all it's information. Usually it is better to use the copy() function for polymorphism. More... | |
~HeatFluidCoupledParticle () override=default | |
HeatFluidCoupledParticle destructor, needs to be implemented and checked if it removes tangential spring information. More... | |
HeatFluidCoupledParticle * | copy () const override |
HeatFluidCoupledParticle copy method. Use copy constructor of this HeatFluidCoupledParticle to create a copy on the heap, useful for polymorphism. More... | |
void | write (std::ostream &os) const override |
HeatFluidCoupledParticle write function, writes HeatFluidCoupledParticle information to the given output-stream, for example a restart-file. More... | |
std::string | getName () const override |
void | read (std::istream &is) override |
HeatFluidCoupledParticle read function, reads in the information for this HeatFluidCoupledParticle from the given input-stream, for example a restart file. More... | |
Mdouble | getLiquidVolume () const |
Returns the volume of the Liquid. More... | |
void | setLiquidVolume (Mdouble liquidVolume) |
Sets the volume of the Liquid. More... | |
void | addLiquidVolume (Mdouble liquidVolume) |
Adds the volume of the Liquid. More... | |
unsigned | getNumberOfFieldsVTK () const override |
std::string | getTypeVTK (unsigned i) const override |
std::string | getNameVTK (unsigned i) const override |
std::vector< Mdouble > | getFieldVTK (unsigned i) const override |
bool | isSphericalParticle () const override |
void | actionsAfterTimeStep () override |
![]() | |
ThermalParticle () | |
Basic Particle constructor, creates a particle at (0,0,0) with radius, mass and inertia equal to 1. More... | |
ThermalParticle (const ThermalParticle &p) | |
Particle copy constructor, which accepts as input a reference to a Particle. It creates a copy of this Particle and all it's information. Usually it is better to use the copy() function for polymorphism. More... | |
~ThermalParticle () override=default | |
Particle destructor, needs to be implemented and checked if it removes tangential spring information. More... | |
ThermalParticle * | copy () const override |
Particle copy method. It calls to copy constructor of this Particle, useful for polymorfism. More... | |
void | write (std::ostream &os) const override |
std::string | getName () const override |
Returns the name of the object. More... | |
void | read (std::istream &is) override |
Mdouble | getTemperature () const |
void | setTemperature (Mdouble temperature) |
void | addTemperature (Mdouble temperature) |
void | setTemperatureDependentDensity (const std::function< double(double)> &temperatureDependentDensity) |
const std::function< double(double)> & | getTemperatureDependentDensity () const |
const std::function< double(double)> & | getTimeDependentTemperature () const |
void | setTimeDependentTemperature (const std::function< double(double)> &timeDependentTemperature) |
void | actionsAfterTimeStep () override |
bool | isSphericalParticle () const override |
![]() | |
BaseParticle () | |
Basic Particle constructor, creates an Particle at (0,0,0) with radius, mass and inertia equal to 1. More... | |
BaseParticle (const BaseParticle &p) | |
Particle copy constructor, which accepts as input a reference to a Particle. It creates a copy of this Particle and all it's information. Usually it is better to use the copy() function for polymorphism. More... | |
BaseParticle (const ParticleSpecies *s) | |
~BaseParticle () override | |
Particle destructor, needs to be implemented and checked if it removes tangential spring information. More... | |
virtual Mdouble | getVolume () const |
Get Particle volume function, which required a reference to the Species vector. It returns the volume of the Particle. More... | |
void | fixParticle () |
Fix Particle function. It fixes a Particle by setting its inverse mass and inertia and velocities to zero. More... | |
bool | isFixed () const override |
Is fixed Particle function. It returns whether a Particle is fixed or not, by checking its inverse Mass. More... | |
bool | isMPIParticle () const |
Indicates if this particle is a ghost in the MPI domain. More... | |
void | setMPIParticle (bool flag) |
Flags the mpi particle status. More... | |
bool | isInMPIDomain () |
Indicates if the particle is in the communication zone of the mpi domain. More... | |
void | setInMPIDomain (bool flag) |
Flags the status of the particle if wether it is in the communication zone or not. More... | |
bool | isInPeriodicDomain () const |
Indicates if the particle is in the periodic boundary communication zone. More... | |
void | setInPeriodicDomain (bool flag) |
Flags the status of the particle wether it is in the periodic communication zone or not. More... | |
bool | isPeriodicGhostParticle () const |
Indicates if this particle is a ghost in the periodic boundary. More... | |
void | setPeriodicGhostParticle (bool flag) |
Flags the status of the particle to be a ghost in periodic boundary or not. More... | |
bool | isMaserParticle () const |
Indicates if this particle belongs to the maser boundary. More... | |
void | setMaserParticle (bool flag) |
Flags the status of the particle if it belongs to the maser boundary or not. More... | |
void | setCommunicationComplexity (unsigned complexity) |
Set the communication complexity of the particle. More... | |
unsigned | getCommunicationComplexity () |
Obtains the communication complexity of the particle. More... | |
void | setPeriodicComplexity (std::vector< int > complexity) |
Set the periodic communication complexity of the particle. More... | |
void | setPeriodicComplexity (int index, int value) |
Set the periodic communication complexity of the particle. More... | |
const std::vector< int > & | getPeriodicComplexity () |
Obtains the periodic communication complexity of the particle. More... | |
void | setPreviousPeriodicComplexity (std::vector< int > complexity) |
Set the previous periodic communication complexity of the paritcle. More... | |
const std::vector< int > & | getPreviousPeriodicComplexity () const |
Sets the previous periodic communication complexity of the particle. More... | |
int | getPeriodicComplexity (int index) |
Gets the periodic communication complexity of a certain boundary. More... | |
void | unfix () |
Unfix Particle function, which required a reference to the Species vector. It unfixes a Particle by computing the Particles mass and inertia. More... | |
virtual void | oldRead (std::istream &is) |
virtual void | setInfo (Mdouble info) |
Sets some user-defined information about this object (by default, species ID). More... | |
virtual Mdouble | getInfo () const |
Returns some user-defined information about this object (by default, species ID). More... | |
void | printHGrid (std::ostream &os) const |
Adds particle's HGrid level and cell coordinates to an ostream. More... | |
unsigned int | getHGridLevel () const |
Returns particle's HGrid level. More... | |
BaseParticle * | getHGridNextObject () const |
Returns pointer to next object in particle's HGrid level & cell. More... | |
BaseParticle * | getHGridPrevObject () const |
Returns pointer to previous object in particle's HGrid level & cell. More... | |
int | getHGridX () const |
Returns particle's HGrid cell X-coordinate. More... | |
int | getHGridY () const |
Returns particle's HGrid cell Y-coordinate. More... | |
int | getHGridZ () const |
Returns particle's HGrid cell Z-coordinate. More... | |
MatrixSymmetric3D | getInvInertia () const |
Returns the inverse of the particle's inertia tensor. More... | |
Mdouble | getInvMass () const override |
Returns the inverse of the particle's mass. More... | |
Mdouble | getCurvature (const Vec3D &labFixedCoordinates) const override |
Mdouble | getKineticEnergy () const |
Calculates the particle's translational kinetic energy. More... | |
Mdouble | getRotationalEnergy () const |
Calculates the particle's rotational kinetic energy. More... | |
Mdouble | getGravitationalEnergy () const |
Calculates the particle's gravitational energy. More... | |
Mdouble | getMass () const |
Returns the particle's mass. More... | |
Mdouble | getSurfaceArea () const |
Vec3D | getMomentum () const |
MatrixSymmetric3D | getInertia () const |
Vec3D | getAngularMomentum () const |
BaseParticle * | getPeriodicFromParticle () const |
Returns the 'original' particle this one's a periodic copy of. More... | |
Mdouble | getRadius () const |
Returns the particle's radius. More... | |
Mdouble | getMaxInteractionRadius () const |
Returns the particle's interaction radius, which might be different from radius_ (e.g., when dealing with wet particles) More... | |
Mdouble | getInteractionDistance (const BaseInteractable *i) const |
Returns the interactionDistance_ of the mixed species of this particle and the particle or wall i. More... | |
Mdouble | getSumOfInteractionRadii (const BaseParticle *particle) const |
returns the sum of the radii plus the interactionDistance More... | |
Mdouble | getWallInteractionRadius (const BaseWall *wall) const |
returns the radius plus the interactionDistance More... | |
const Vec3D & | getDisplacement () const |
Returns the particle's displacement relative to the previous time step. More... | |
const Vec3D & | getPreviousPosition () const |
Returns the particle's position in the previous time step. More... | |
const Vec3D | getDisplacement2 (Mdouble xmin, Mdouble xmax, Mdouble ymin, Mdouble ymax, Mdouble zmin, Mdouble zmax, Mdouble t) const |
virtual void | setInertia () |
void | setInertia (MatrixSymmetric3D inertia) |
Sets the particle's inertia_ (and adjusts invInertia_ accordingly) More... | |
void | setInverseInertia (MatrixSymmetric3D inverseInertia) |
Sets the particle's inertia_ (and adjusts invInertia_ accordingly) More... | |
void | setInfiniteInertia () |
Sets the particle's inertia_ to 'infinite' (1e20) and its invInertia_ to 0. More... | |
void | setPeriodicFromParticle (BaseParticle *p) |
Assigns the pointer to the 'original' particle this one's a periodic copy of (used in periodic boundary condition implementations). More... | |
void | setHGridX (const int x) |
Sets the particle's HGrid cell X-coordinate. More... | |
void | setHGridY (const int y) |
Sets the particle's HGrid cell Y-coordinate. More... | |
void | setHGridZ (const int z) |
Sets the particle's HGrid cell Z-coordinate. More... | |
void | setHGridLevel (const unsigned int level) |
Sets the particle's HGrid level. More... | |
void | setHGridNextObject (BaseParticle *p) |
Sets the pointer to the next object in the particle's HGrid cell & level. More... | |
void | setHGridPrevObject (BaseParticle *p) |
Sets the pointer to the previous object in the particle's HGrid cell & level. More... | |
virtual void | setRadius (Mdouble radius) |
Sets the particle's radius_ (and adjusts the mass_ accordingly, based on the particle's species) More... | |
virtual Vec3D | getAxes () const |
Only ustilised in case of superquadric particles. Had to create a virtual function to allow function access in writeVTK function in the particle handler. More... | |
virtual Mdouble | getExponentEps1 () const |
Only ustilised in case of superquadric particles. Had to create a virtual function to allow function access in writeVTK function in the particle handler. More... | |
virtual Mdouble | getExponentEps2 () const |
Only ustilised in case of superquadric particles. Had to create a virtual function to allow function access in writeVTK function in the particle handler. More... | |
virtual void | setAxes (const Vec3D &axes) |
Only ustilised in case of superquadric particles. More... | |
virtual void | setExponents (const Mdouble &eps1, const Mdouble &eps2) |
Only ustilised in case of superquadric particles. More... | |
MERCURY_DEPRECATED void | setMass (Mdouble mass) |
Sets the particle's mass. More... | |
void | setMassForP3Statistics (Mdouble mass) |
Sets the particle's mass This function should not be used, but is necessary to extend the CG toolbox to non-spherical particles. More... | |
void | setDisplacement (const Vec3D &disp) |
Sets the particle's displacement (= difference between current position and that of the previous time step) More... | |
void | setPreviousPosition (const Vec3D &pos) |
Sets the particle's position in the previous time step. More... | |
void | movePrevious (const Vec3D &posMove) |
Adds a vector to the particle's previousPosition_. More... | |
void | accelerate (const Vec3D &vel) |
Increases the particle's velocity_ by the given vector. More... | |
void | angularAccelerate (const Vec3D &angVel) |
Increases the particle's angularVelocity_ by the given vector. More... | |
void | addDisplacement (const Vec3D &addDisp) |
Adds a vector to the particle's displacement_. More... | |
void | setHandler (ParticleHandler *handler) |
Sets the pointer to the particle's ParticleHandler. More... | |
ParticleHandler * | getHandler () const |
Returns pointer to the particle's ParticleHandler. More... | |
BaseInteraction * | getInteractionWith (BaseParticle *P, unsigned timeStamp, InteractionHandler *interactionHandler) override |
Checks if particle is in interaction with given particle P, and if so, returns vector of pointer to the associated BaseInteraction object (else returns empty vector). More... | |
virtual bool | isInContactWith (const BaseParticle *P) const |
Get whether or not this particle is in contact with the given particle. More... | |
void | integrateBeforeForceComputation (double time, double timeStep) |
First step of Velocity Verlet integration. More... | |
void | integrateAfterForceComputation (double time, double timeStep) |
Second step of Velocity Verlet integration. More... | |
unsigned int | getParticleDimensions () const |
Returns the particle's dimensions (either 2 or 3). More... | |
MERCURY_DEPRECATED void | setIndSpecies (unsigned int indSpecies) override |
void | setSpecies (const ParticleSpecies *species) |
const HGridCell & | getHGridCell () const |
virtual void | computeMass (const ParticleSpecies &s) |
Computes the particle's (inverse) mass and inertia. More... | |
![]() | |
BaseInteractable () | |
Default BaseInteractable constructor. More... | |
BaseInteractable (const BaseInteractable &p) | |
Copy constructor. More... | |
~BaseInteractable () override | |
Destructor, it simply destructs the BaseInteractable and all the objects it contains. More... | |
unsigned int | getIndSpecies () const |
Returns the index of the species associated with the interactable object. More... | |
const ParticleSpecies * | getSpecies () const |
Returns a pointer to the species of this BaseInteractable. More... | |
void | setSpecies (const ParticleSpecies *species) |
Sets the species of this BaseInteractable. More... | |
const Vec3D & | getForce () const |
Returns the force on this BaseInteractable. More... | |
const Vec3D & | getTorque () const |
Returns the torque on this BaseInteractable. More... | |
void | setForce (const Vec3D &force) |
Sets the force on this BaseInteractable. More... | |
void | setTorque (const Vec3D &torque) |
Sets the torque on this BaseInteractable. More... | |
void | addForce (const Vec3D &addForce) |
Adds an amount to the force on this BaseInteractable. More... | |
void | addTorque (const Vec3D &addTorque) |
Adds an amount to the torque on this BaseInteractable. More... | |
void | resetForceTorque (int numberOfOMPthreads) |
void | sumForceTorqueOMP () |
const Vec3D & | getPosition () const |
Returns the position of this BaseInteractable. More... | |
const Quaternion & | getOrientation () const |
Returns the orientation of this BaseInteractable. More... | |
void | setPosition (const Vec3D &position) |
Sets the position of this BaseInteractable. More... | |
void | setOrientationViaNormal (Vec3D normal) |
Sets the orientation of this BaseInteractable by defining the vector that results from the rotation of the (1,0,0) vector. More... | |
void | setOrientationViaEuler (Vec3D eulerAngle) |
Sets the orientation of this BaseInteractable by defining the euler angles. More... | |
void | setOrientation (const Quaternion &orientation) |
Sets the orientation of this BaseInteractable. More... | |
virtual void | move (const Vec3D &move) |
Moves this BaseInteractable by adding an amount to the position. More... | |
virtual void | rotate (const Vec3D &angularVelocityDt) |
Rotates this BaseInteractable. More... | |
const std::vector < BaseInteraction * > & | getInteractions () const |
Returns a list of interactions which belong to this interactable. More... | |
void | addInteraction (BaseInteraction *I) |
Adds an interaction to this BaseInteractable. More... | |
bool | removeInteraction (BaseInteraction *I) |
Removes an interaction from this BaseInteractable. More... | |
void | copyInteractionsForPeriodicParticles (const BaseInteractable &p) |
Copies interactions to this BaseInteractable whenever a periodic copy made. More... | |
void | setVelocity (const Vec3D &velocity) |
set the velocity of the BaseInteractable. More... | |
void | setAngularVelocity (const Vec3D &angularVelocity) |
set the angular velocity of the BaseInteractble. More... | |
void | addVelocity (const Vec3D &velocity) |
adds an increment to the velocity. More... | |
void | addAngularVelocity (const Vec3D &angularVelocity) |
add an increment to the angular velocity. More... | |
virtual const Vec3D & | getVelocity () const |
Returns the velocity of this interactable. More... | |
virtual const Vec3D & | getAngularVelocity () const |
Returns the angular velocity of this interactable. More... | |
void | setPrescribedPosition (const std::function< Vec3D(double)> &prescribedPosition) |
Allows the position of an infinite mass interactable to be prescribed. More... | |
void | applyPrescribedPosition (double time) |
Computes the position from the user defined prescribed position function. More... | |
void | setPrescribedVelocity (const std::function< Vec3D(double)> &prescribedVelocity) |
Allows the velocity of an infinite mass interactable to be prescribed. More... | |
void | applyPrescribedVelocity (double time) |
Computes the velocity from the user defined prescribed velocity function. More... | |
void | setPrescribedOrientation (const std::function< Quaternion(double)> &prescribedOrientation) |
Allows the orientation of the infinite mass interactbale to be prescribed. More... | |
void | applyPrescribedOrientation (double time) |
Computes the orientation from the user defined prescribed orientation function. More... | |
void | setPrescribedAngularVelocity (const std::function< Vec3D(double)> &prescribedAngularVelocity) |
Allows the angular velocity of the infinite mass interactable to be prescribed. More... | |
void | applyPrescribedAngularVelocity (double time) |
Computes the angular velocity from the user defined prescribed angular velocity. More... | |
virtual const Vec3D | getVelocityAtContact (const Vec3D &contact) const |
Returns the velocity at the contact point, use by many force laws. More... | |
void | integrateBeforeForceComputation (double time, double timeStep) |
This is part of integrate routine for objects with infinite mass. More... | |
void | integrateAfterForceComputation (double time, double timeStep) |
This is part of the integration routine for objects with infinite mass. More... | |
![]() | |
BaseObject ()=default | |
Default constructor. More... | |
BaseObject (const BaseObject &p)=default | |
Copy constructor, copies all the objects BaseObject contains. More... | |
virtual | ~BaseObject ()=default |
virtual destructor More... | |
virtual void | moveInHandler (unsigned int index) |
Except that it is virtual, it does the same thing as setIndex() does. More... | |
void | setIndex (unsigned int index) |
Allows one to assign an index to an object in the handler/container. More... | |
void | setId (unsigned long id) |
Assigns a unique identifier to each object in the handler (container) which remains constant even after the object is deleted from the container/handler. More... | |
unsigned int | getIndex () const |
Returns the index of the object in the handler. More... | |
unsigned int | getId () const |
Returns the unique identifier of any particular object. More... | |
void | setGroupId (unsigned groupId) |
unsigned | getGroupId () const |
Private Member Functions | |
double | f1 (double liquidVolume, double temperature) |
f1 is used in Runge–Kutta method. More... | |
double | f2 (double liquidVolume, double temperature) |
f2 is used in Runge–Kutta method. More... | |
Private Attributes | |
Mdouble | liquidVolume_ |
Additional Inherited Members | |
![]() | |
Mdouble | temperature_ |
![]() | |
Mdouble | radius_ |
Mdouble | invMass_ |
Particle radius_. More... | |
MatrixSymmetric3D | invInertia_ |
Inverse Particle mass (for computation optimization) More... | |
Class that implements particles which store both temperature/heat capacity and liquid content which is adapted for the CFD-DEM studies.
In some CFD-DEM studies, a drying process needs to be simulated. Therefore, we selected the drying model of Azmir et al. (2018), which considers heat and mass transfer. Their model considers convection, conduction, radiation, and evaporation. Conduction and convection had already been implemented in ThermalParticle class. This class implements evaporation model.
Definition at line 41 of file HeatFluidCoupledParticle.h.
|
inline |
HeatFluidCoupledParticle constructor creates a HeatFluidCoupledParticle at (0,0,0) with radius, mass and inertia equal to 1.
Definition at line 48 of file HeatFluidCoupledParticle.h.
References liquidVolume_.
Referenced by copy().
|
inline |
HeatFluidCoupledParticle copy constructor, which accepts as input a reference to a HeatFluidCoupledParticle. It creates a copy of this HeatFluidCoupledParticle and all it's information. Usually it is better to use the copy() function for polymorphism.
Constructor that copies most of the properties of the given particle. Please note that not everything is copied, for example the position in the HGrid is not determined yet by the end of this constructor. It also does not copy the interactions and the pointer to the handler that handles this particle. Use with care.
[in,out] | p | Reference to the HeatFluidCoupledParticle this one should become a copy of. |
Definition at line 64 of file HeatFluidCoupledParticle.h.
References liquidVolume_.
|
overridedefault |
HeatFluidCoupledParticle destructor, needs to be implemented and checked if it removes tangential spring information.
Destructor that asks the ParticleHandler to check if this was the smallest or largest particle and adjust itself accordingly.
|
overridevirtual |
Reimplemented from BaseParticle.
Definition at line 85 of file HeatFluidCoupledParticle.cc.
References f1(), f2(), BaseHandler< T >::getDPMBase(), BaseParticle::getHandler(), DPMBase::getTimeStep(), liquidVolume_, and ThermalParticle::temperature_.
|
inline |
Adds the volume of the Liquid.
Definition at line 134 of file HeatFluidCoupledParticle.h.
References liquidVolume_.
|
inlineoverridevirtual |
HeatFluidCoupledParticle copy method. Use copy constructor of this HeatFluidCoupledParticle to create a copy on the heap, useful for polymorphism.
Implements BaseParticle.
Definition at line 81 of file HeatFluidCoupledParticle.h.
References HeatFluidCoupledParticle().
|
private |
f1 is used in Runge–Kutta method.
Definition at line 102 of file HeatFluidCoupledParticle.cc.
References mathsFunc::exp(), BaseParticle::getMass(), BaseInteractable::getSpecies(), BaseParticle::getSurfaceArea(), mathsFunc::log(), and constants::R.
Referenced by actionsAfterTimeStep().
|
private |
f2 is used in Runge–Kutta method.
Definition at line 114 of file HeatFluidCoupledParticle.cc.
References mathsFunc::exp(), BaseParticle::getMass(), BaseInteractable::getSpecies(), BaseParticle::getSurfaceArea(), mathsFunc::log(), and constants::R.
Referenced by actionsAfterTimeStep().
Reimplemented from BaseParticle.
Definition at line 65 of file HeatFluidCoupledParticle.cc.
References BaseInteractable::getInteractions(), LiquidMigrationWilletInteraction::getLiquidBridgeVolume(), liquidVolume_, and ThermalParticle::temperature_.
|
inline |
Returns the volume of the Liquid.
Definition at line 118 of file HeatFluidCoupledParticle.h.
References liquidVolume_.
|
inlineoverridevirtual |
Returns the name of the object; in this case "HeatFluidCoupledParticle".
Reimplemented from BaseParticle.
Definition at line 103 of file HeatFluidCoupledParticle.h.
|
overridevirtual |
Reimplemented from BaseParticle.
Definition at line 53 of file HeatFluidCoupledParticle.cc.
|
inlineoverridevirtual |
Reimplemented from BaseParticle.
Definition at line 139 of file HeatFluidCoupledParticle.h.
|
inlineoverridevirtual |
Reimplemented from BaseParticle.
Definition at line 144 of file HeatFluidCoupledParticle.h.
|
inlineoverridevirtual |
Reimplemented from BaseParticle.
Definition at line 153 of file HeatFluidCoupledParticle.h.
|
overridevirtual |
HeatFluidCoupledParticle read function, reads in the information for this HeatFluidCoupledParticle from the given input-stream, for example a restart file.
Particle read function. Has an std::istream as argument, from which it extracts the radius_, invMass_ and invInertia_, respectively. From these the mass_ and inertia_ are deduced. An additional set of properties is read through the call to the parent's method ThermalParticle::read().
[in,out] | is | input stream with particle properties. |
Reimplemented from BaseParticle.
Definition at line 46 of file HeatFluidCoupledParticle.cc.
References liquidVolume_, and ThermalParticle::read().
|
inline |
Sets the volume of the Liquid.
Definition at line 126 of file HeatFluidCoupledParticle.h.
References liquidVolume_.
|
inlineoverridevirtual |
HeatFluidCoupledParticle write function, writes HeatFluidCoupledParticle information to the given output-stream, for example a restart-file.
HeatFluidCoupledParticle print method, which accepts an os std::ostream as input. It prints human readable HeatFluidCoupledParticle information to the std::ostream.
[in,out] | os | stream to which the info is written. |
Reimplemented from BaseParticle.
Definition at line 93 of file HeatFluidCoupledParticle.h.
References liquidVolume_, and ThermalParticle::write().
|
private |
Definition at line 166 of file HeatFluidCoupledParticle.h.
Referenced by actionsAfterTimeStep(), addLiquidVolume(), getFieldVTK(), getLiquidVolume(), HeatFluidCoupledParticle(), read(), setLiquidVolume(), and write().