ECOGEN 4.0
Evolutive, Compressible, Open, Genuine, Easy, N-phase
Loading...
Searching...
No Matches
Cell Class Reference

Base class for a mesh cell. More...

#include <Cell.h>

Inherited by CellGhost, and CellO2.

Public Member Functions

 Cell ()
 Basic Cell constructor for a non AMR cell.
 
 Cell (int lvl)
 Cell constructor for an AMR cell.
 
virtual ~Cell ()
 
void addCellInterface (CellInterface *cellInterface)
 Add a cell interface to current cell.
 
void deleteCellInterface (CellInterface *cellInterface)
 Delete the cell interface of current cell.
 
void associateExtVar (Model *mod, Gradient *grad)
 Associate external variables.
 
virtual void allocate (const std::vector< AddPhys * > &addPhys)
 Memory allocation of cell attributes.
 
void allocateEos ()
 
void fill (std::vector< GeometricalDomain * > &domains, const int &)
 Filling cell properties using a physical domain.
 
virtual void copyPhase (const int &phaseNumber, Phase *phase)
 
void copyMixture (Mixture *mixture)
 
void setToZeroCons ()
 
void setToZeroConsGlobal ()
 
void timeEvolution (const double &dt, Symmetry *symmetry)
 
void timeEvolutionAddPhys (const double &dt)
 
void buildPrim ()
 
void buildCons ()
 
void correctionEnergy ()
 
void sourceTermIntegration (const double &)
 
void printPhasesMixture (std::ofstream &fileStream) const
 
virtual void completeFulfillState ()
 
virtual void fulfillState (Prim=vecPhases)
 
void fulfillStateRestart ()
 
void localProjection (const Coord &normal, const Coord &tangent, const Coord &binormal, Prim=vecPhases)
 
virtual void reverseProjection (const Coord &normal, const Coord &tangent, const Coord &binormal)
 
virtual void copyInCell (Cell &) const
 
void copyVec (Phase **vecPhases, Mixture *mixture, Transport *vecTransports)
 
void deleteInterface (const int &b)
 
void prepareAddPhys ()
 
Coord selectVector (Variable nameVector, int num=0, int subscript=-1) const
 
void setVector (Variable nameVector, const Coord &value, int num=0, int subscript=-1)
 
const QuantitiesAddPhysgetQPA (int &numQPA) const
 
const CoordgetGradTk (int &numPhase, int &numAddPhys) const
 Allow to recover an additional physical quantity.
 
void setGradTk (int &numPhase, int &numAddPhys, double *buffer, int &counter)
 
void addNonConsAddPhys (AddPhys &addPhys, Symmetry *symmetry)
 
void reinitializeColorFunction (const int &numTransport, const int &numPhase)
 Re-initialize the color function (transport) with alpha.
 
void computeGradients (std::vector< Coord > &grads, std::vector< Variable > &nameVariables, std::vector< int > &numPhases)
 Compute gradients (temperature, velocity, density) of a cell.
 
int getCellInterfacesSize () const
 
CellInterfacegetCellInterface (const int &b)
 
void setCellInterface (const int &b, CellInterface *cellInterface)
 
virtual PhasegetPhase (const int &phaseNumber, Prim=vecPhases) const
 
virtual Phase ** getPhases (Prim=vecPhases) const
 
virtual MixturegetMixture (Prim=vecPhases) const
 
FluxgetCons () const
 
void setCons (Flux *cons)
 
const CoordgetPosition () const
 
const CoordgetSize () const
 
const double & getSizeX () const
 
const double & getSizeY () const
 
const double & getSizeZ () const
 
void setElement (Element *element, const int &numCell)
 
ElementgetElement () const
 
virtual void setTransport (double value, int &numTransport, Prim=vecPhases)
 
virtual TransportgetTransport (const int &numTransport, Prim=vecPhases) const
 
virtual TransportgetTransports (Prim=vecPhases) const
 
TransportgetConsTransport (const int &numTransport) const
 
void setConsTransport (double value, const int &numTransport)
 
std::vector< QuantitiesAddPhys * > & getVecQuantitiesAddPhys ()
 
const int & getNumberPhases () const
 
const int & getNumberTransports () const
 
double getDensityGradient ()
 
ModelgetModel ()
 
CoordgetVelocity ()
 
const CoordgetVelocity () const
 
void setWall (bool wall)
 
bool getWall () const
 
double selectScalar (Variable nameVariable, int num=0) const
 Select a specific scalar variable.
 
virtual void allocateSecondOrderBuffersAndGradientVectors (Phase **, Mixture *)
 Compute global variable buffers (min, max, etc.) and initialize speficic gradient vectors for 2nd-order scheme on unstructured mesh.
 
virtual void computeGradientsO2 ()
 Compute gradients for 2nd-order scheme on unstructured mesh.
 
virtual void limitGradientsO2 (Limiter &)
 
virtual void computeLocalSlopes (CellInterface &)
 Compute slopes for 2nd-order scheme on unstructured mesh.
 
virtual void computeLocalSlopes (CellInterface &, Limiter &, Limiter &, Limiter &, Limiter &, double &, double &, double &, double &)
 
virtual void computeLocalSlopesLimite (CellInterface &, Limiter &, Limiter &, Limiter &, Limiter &, double &)
 
virtual PhasegetSlopes () const
 
virtual TransportgetSlopesTransport () const
 
virtual void saveCons ()
 
virtual void getBackCons ()
 
virtual void predictionOrdre2 (const double &, Symmetry *)
 
double distance (Cell *c)
 
double distanceX (Cell *c)
 
double distanceY (Cell *c)
 
double distanceZ (Cell *c)
 
double distance (CellInterface *b)
 
double distanceX (CellInterface *b)
 
double distanceY (CellInterface *b)
 
double distanceZ (CellInterface *b)
 
bool traverseObjet (const GeometricObject &objet) const
 
void printInfo () const
 
double getPsat ()
 Compute saturation pressure for a liquid/vapor fluid mixture (first phase is considered predominant -> not generalized be careful)
 
bool printGnuplotAMR (std::ofstream &fileStream, const int &dim, GeometricObject *objet=0, bool recordPsat=false)
 
void computeVolumePhaseK (double &integration, const int &numPhase)
 
void computeMass (double &mass, double &alphaRef)
 
void computeTotalMass (double &mass)
 
void computeTotalEnergy (double &totalEnergy)
 
void lookForPmax (double *pMax, double *pMaxWall)
 
void setToZeroXi ()
 
void setToZeroConsXi ()
 
void timeEvolutionXi ()
 
void chooseRefine (const double &xiSplit, const int &nbCellsY, const int &nbCellsZ, const std::vector< AddPhys * > &addPhys, int &nbCellsTotalAMR)
 
void chooseUnrefine (const double &xiJoin, int &nbCellsTotalAMR)
 
void refineCellAndCellInterfaces (const int &nbCellsY, const int &nbCellsZ, const std::vector< AddPhys * > &addPhys, const bool &refineExternalCellInterfaces)
 
virtual void createChildCell (const int &lvl)
 
void unrefineCellAndCellInterfaces ()
 
void averageChildrenInParent ()
 
bool lvlNeighborTooHigh ()
 
bool lvlNeighborTooLow ()
 
void buildLvlCellsAndLvlInternalCellInterfacesArrays (std::vector< Cell * > *cellsLvl, std::vector< CellInterface * > *cellInterfacesLvl)
 
const int & getLvl () const
 
const bool & getSplit () const
 
const double & getXi () const
 
void setXi (double value)
 
void addFluxXi (double value)
 
int getNumberCellsChildren ()
 
CellgetCellChild (const int &num)
 
std::vector< Cell * > * getChildVector ()
 
virtual void pushBackSlope ()
 
virtual void popBackSlope ()
 
virtual int getRankOfNeighborCPU () const
 
virtual void setRankOfNeighborCPU (int)
 
void fillBufferPrimitives (double *buffer, int &counter, const int &lvl, const int &neighbour, Prim type=vecPhases) const
 
void getBufferPrimitives (double *buffer, int &counter, const int &lvl, Eos **eos, Prim type=vecPhases)
 
void fillBufferVector (double *buffer, int &counter, const int &lvl, const int &neighbour, const int &dim, Variable nameVector, int num=0, int index=-1) const
 
void getBufferVector (double *buffer, int &counter, const int &lvl, const int &dim, Variable nameVector, int num=0, int index=-1)
 
void fillBufferTransports (double *buffer, int &counter, const int &lvl, const int &neighbour) const
 
void getBufferTransports (double *buffer, int &counter, const int &lvl)
 
virtual void fillBufferSlopes (double *, int &, const int &, const int &) const
 
virtual void getBufferSlopes (double *, int &, const int &)
 
virtual bool isCellGhost () const
 
bool hasNeighboringGhostCellOfCPUneighbour (const int &neighbour) const
 
int numberOfNeighboringGhostCellsOfCPUneighbour (const int &neighbour) const
 
virtual GradPhasegetGradPhase (const int &) const
 
virtual GradMixturegetGradMixture () const
 
virtual GradTransportgetGradTransport (const int &) const
 
void chooseRefineDeraffineGhost (const int &nbCellsY, const int &nbCellsZ, const std::vector< AddPhys * > &addPhys, std::vector< Cell * > *cellsLvlGhost)
 
void refineCellAndCellInterfacesGhost (const int &nbCellsY, const int &nbCellsZ, const std::vector< AddPhys * > &addPhys)
 
void unrefineCellAndCellInterfacesGhost ()
 
void fillBufferXi (double *buffer, int &counter, const int &lvl, const int &neighbour) const
 
void getBufferXi (double *buffer, int &counter, const int &lvl)
 
void fillBufferSplit (bool *buffer, int &counter, const int &lvl, const int &neighbour) const
 
void getBufferSplit (bool *buffer, int &counter, const int &lvl)
 
void fillNumberElementsToSendToNeighbour (int &numberElementsToSendToNeighbor, int &numberSlopesToSendToNeighbor, const int &lvl, const int &neighbour, int numberNeighboursOfCPUneighbour)
 
void fillDataToSend (std::vector< double > &dataToSend, std::vector< int > &dataSplitToSend, const int &lvl) const
 
void getDataToReceiveAndRefine (std::vector< double > &dataToReceive, std::vector< int > &dataSplitToReceive, const int &lvl, Eos **eos, int &counter, int &counterSplit, const int &nbCellsY, const int &nbCellsZ, const std::vector< AddPhys * > &addPhys)
 
void computeLoad (double &load, int lvl) const
 
void computeLvlMax (int &lvlMax) const
 
void clearExternalCellInterfaces (const int &nbCellsY, const int &nbCellsZ)
 
void updatePointersInternalCellInterfaces ()
 
void updateNbCellsTotalAMR (int &nbCellsTotalAMR)
 

Protected Attributes

bool m_wall
 
Phase ** m_vecPhases
 
Mixturem_mixture
 
Transportm_vecTransports
 
Fluxm_cons
 
Transportm_consTransports
 
Elementm_element
 
std::vector< CellInterface * > m_cellInterfaces
 
std::vector< QuantitiesAddPhys * > m_vecQuantitiesAddPhys
 
int m_lvl
 
double m_xi
 
double m_consXi
 
bool m_split
 
std::vector< Cell * > m_childrenCells
 
std::vector< CellInterface * > m_childrenInternalCellInterfaces
 

Detailed Description

Base class for a mesh cell.

Constructor & Destructor Documentation

◆ Cell() [1/2]

Cell::Cell ( )

Basic Cell constructor for a non AMR cell.

◆ Cell() [2/2]

Cell::Cell ( int  lvl)

Cell constructor for an AMR cell.

Parameters
lvllevel of current AMR cell

◆ ~Cell()

Cell::~Cell ( )
virtual

Member Function Documentation

◆ addCellInterface()

void Cell::addCellInterface ( CellInterface cellInterface)

Add a cell interface to current cell.

Parameters
cellInterfacepointer to added cell interface

◆ addFluxXi()

void Cell::addFluxXi ( double  value)

Add xi cell flux

◆ addNonConsAddPhys()

void Cell::addNonConsAddPhys ( AddPhys addPhys,
Symmetry symmetry 
)

◆ allocate()

void Cell::allocate ( const std::vector< AddPhys * > &  addPhys)
virtual

Memory allocation of cell attributes.

Parameters
addPhysvector of additional physics

Reimplemented in CellO2, CellO2GhostCartesian, and CellO2NS.

◆ allocateEos()

void Cell::allocateEos ( )

◆ allocateSecondOrderBuffersAndGradientVectors()

virtual void Cell::allocateSecondOrderBuffersAndGradientVectors ( Phase **  ,
Mixture  
)
inlinevirtual

Compute global variable buffers (min, max, etc.) and initialize speficic gradient vectors for 2nd-order scheme on unstructured mesh.

Reimplemented in CellO2, and CellO2NS.

◆ associateExtVar()

void Cell::associateExtVar ( Model mod,
Gradient grad 
)

Associate external variables.

Parameters
modpointer to model
gradpointer to gradient method

◆ averageChildrenInParent()

void Cell::averageChildrenInParent ( )

Average variables of children cells in the parent cell, needed for computation of xi.

◆ buildCons()

void Cell::buildCons ( )

◆ buildLvlCellsAndLvlInternalCellInterfacesArrays()

void Cell::buildLvlCellsAndLvlInternalCellInterfacesArrays ( std::vector< Cell * > *  cellsLvl,
std::vector< CellInterface * > *  cellInterfacesLvl 
)

Build new arrays of cells and cell interfaces for level (lvl+1), only internal cell interfaces are added here

◆ buildPrim()

void Cell::buildPrim ( )

◆ chooseRefine()

void Cell::chooseRefine ( const double &  xiSplit,
const int &  nbCellsY,
const int &  nbCellsZ,
const std::vector< AddPhys * > &  addPhys,
int &  nbCellsTotalAMR 
)

Choice for refinement of parent cell

◆ chooseRefineDeraffineGhost()

void Cell::chooseRefineDeraffineGhost ( const int &  nbCellsY,
const int &  nbCellsZ,
const std::vector< AddPhys * > &  addPhys,
std::vector< Cell * > *  cellsLvlGhost 
)

Choice for refinement, unrefinement of the ghost parent cell + Update of ghost cell vector for lvl+1

◆ chooseUnrefine()

void Cell::chooseUnrefine ( const double &  xiJoin,
int &  nbCellsTotalAMR 
)

Choice for unrefinement of parent cell

◆ clearExternalCellInterfaces()

void Cell::clearExternalCellInterfaces ( const int &  nbCellsY,
const int &  nbCellsZ 
)

◆ completeFulfillState()

void Cell::completeFulfillState ( )
virtual

◆ computeGradients()

void Cell::computeGradients ( std::vector< Coord > &  grads,
std::vector< Variable > &  nameVariables,
std::vector< int > &  numPhases 
)

Compute gradients (temperature, velocity, density) of a cell.

Parameters
gradsArray of desired gradients, e.g. for temperature each component represents phase temperature and for velocity each component represents the gradient of a velocity component (grad(u), grad(v), grad(w))
nameVariablesName of the variable for which the gradient is calculated
numPhasesPhases number's

◆ computeGradientsO2()

virtual void Cell::computeGradientsO2 ( )
inlinevirtual

Compute gradients for 2nd-order scheme on unstructured mesh.

Reimplemented in CellO2, and CellO2NS.

◆ computeLoad()

void Cell::computeLoad ( double &  load,
int  lvl 
) const

◆ computeLocalSlopes() [1/2]

virtual void Cell::computeLocalSlopes ( CellInterface )
inlinevirtual

Compute slopes for 2nd-order scheme on unstructured mesh.

Reimplemented in CellO2, and CellO2NS.

◆ computeLocalSlopes() [2/2]

virtual void Cell::computeLocalSlopes ( CellInterface ,
Limiter ,
Limiter ,
Limiter ,
Limiter ,
double &  ,
double &  ,
double &  ,
double &   
)
inlinevirtual

Reimplemented in CellO2, CellO2Cartesian, and CellO2GhostCartesian.

◆ computeLocalSlopesLimite()

virtual void Cell::computeLocalSlopesLimite ( CellInterface ,
Limiter ,
Limiter ,
Limiter ,
Limiter ,
double &   
)
inlinevirtual

Does nothing for first order cells

Reimplemented in CellO2, and CellO2Cartesian.

◆ computeLvlMax()

void Cell::computeLvlMax ( int &  lvlMax) const

◆ computeMass()

void Cell::computeMass ( double &  mass,
double &  alphaRef 
)

◆ computeTotalEnergy()

void Cell::computeTotalEnergy ( double &  totalEnergy)

◆ computeTotalMass()

void Cell::computeTotalMass ( double &  mass)

◆ computeVolumePhaseK()

void Cell::computeVolumePhaseK ( double &  integration,
const int &  numPhase 
)

◆ copyInCell()

virtual void Cell::copyInCell ( Cell ) const
inlinevirtual

◆ copyMixture()

void Cell::copyMixture ( Mixture mixture)

◆ copyPhase()

void Cell::copyPhase ( const int &  phaseNumber,
Phase phase 
)
virtual

Reimplemented in CellO2.

◆ copyVec()

void Cell::copyVec ( Phase **  vecPhases,
Mixture mixture,
Transport vecTransports 
)

◆ correctionEnergy()

void Cell::correctionEnergy ( )

◆ createChildCell()

void Cell::createChildCell ( const int &  lvl)
virtual

Create a child cell (not initialized)

Reimplemented in CellO2, CellGhost, CellO2Cartesian, and CellO2GhostCartesian.

◆ deleteCellInterface()

void Cell::deleteCellInterface ( CellInterface cellInterface)

Delete the cell interface of current cell.

Parameters
cellInterfacepointer to deleted cell interface

◆ deleteInterface()

void Cell::deleteInterface ( const int &  b)

◆ distance() [1/2]

double Cell::distance ( Cell c)

Does nothing for first order cells Distance totale

◆ distance() [2/2]

double Cell::distance ( CellInterface b)

Distance totale

◆ distanceX() [1/2]

double Cell::distanceX ( Cell c)

Distance selon x

◆ distanceX() [2/2]

double Cell::distanceX ( CellInterface b)

Distance selon x

◆ distanceY() [1/2]

double Cell::distanceY ( Cell c)

Distance selon y

◆ distanceY() [2/2]

double Cell::distanceY ( CellInterface b)

Distance selon y

◆ distanceZ() [1/2]

double Cell::distanceZ ( Cell c)

Distance selon z

◆ distanceZ() [2/2]

double Cell::distanceZ ( CellInterface b)

Distance selon z

◆ fill()

void Cell::fill ( std::vector< GeometricalDomain * > &  domains,
const int &   
)

Filling cell properties using a physical domain.

Parameters
domainsdomain used for filling

◆ fillBufferPrimitives()

void Cell::fillBufferPrimitives ( double *  buffer,
int &  counter,
const int &  lvl,
const int &  neighbour,
Prim  type = vecPhases 
) const

Does nothing for non-ghost cells

◆ fillBufferSlopes()

virtual void Cell::fillBufferSlopes ( double *  ,
int &  ,
const int &  ,
const int &   
) const
inlinevirtual

Reimplemented in CellO2, CellO2NS, and CellO2Cartesian.

◆ fillBufferSplit()

void Cell::fillBufferSplit ( bool *  buffer,
int &  counter,
const int &  lvl,
const int &  neighbour 
) const

◆ fillBufferTransports()

void Cell::fillBufferTransports ( double *  buffer,
int &  counter,
const int &  lvl,
const int &  neighbour 
) const

◆ fillBufferVector()

void Cell::fillBufferVector ( double *  buffer,
int &  counter,
const int &  lvl,
const int &  neighbour,
const int &  dim,
Variable  nameVector,
int  num = 0,
int  index = -1 
) const

◆ fillBufferXi()

void Cell::fillBufferXi ( double *  buffer,
int &  counter,
const int &  lvl,
const int &  neighbour 
) const

◆ fillDataToSend()

void Cell::fillDataToSend ( std::vector< double > &  dataToSend,
std::vector< int > &  dataSplitToSend,
const int &  lvl 
) const

◆ fillNumberElementsToSendToNeighbour()

void Cell::fillNumberElementsToSendToNeighbour ( int &  numberElementsToSendToNeighbor,
int &  numberSlopesToSendToNeighbor,
const int &  lvl,
const int &  neighbour,
int  numberNeighboursOfCPUneighbour 
)

◆ fulfillState()

void Cell::fulfillState ( Prim  = vecPhases)
virtual

Reimplemented in CellO2.

◆ fulfillStateRestart()

void Cell::fulfillStateRestart ( )

◆ getBackCons()

virtual void Cell::getBackCons ( )
inlinevirtual

Does nothing for first order cells

Reimplemented in CellO2.

◆ getBufferPrimitives()

void Cell::getBufferPrimitives ( double *  buffer,
int &  counter,
const int &  lvl,
Eos **  eos,
Prim  type = vecPhases 
)

◆ getBufferSlopes()

virtual void Cell::getBufferSlopes ( double *  ,
int &  ,
const int &   
)
inlinevirtual

Does nothing for first order cells

Reimplemented in CellO2, CellO2Cartesian, CellO2NS, CellO2GhostNS, and CellO2GhostCartesian.

◆ getBufferSplit()

void Cell::getBufferSplit ( bool *  buffer,
int &  counter,
const int &  lvl 
)

◆ getBufferTransports()

void Cell::getBufferTransports ( double *  buffer,
int &  counter,
const int &  lvl 
)

◆ getBufferVector()

void Cell::getBufferVector ( double *  buffer,
int &  counter,
const int &  lvl,
const int &  dim,
Variable  nameVector,
int  num = 0,
int  index = -1 
)

◆ getBufferXi()

void Cell::getBufferXi ( double *  buffer,
int &  counter,
const int &  lvl 
)

◆ getCellChild()

Cell * Cell::getCellChild ( const int &  num)

Return pointer to the corresponding child cell number

◆ getCellInterface()

CellInterface * Cell::getCellInterface ( const int &  b)

◆ getCellInterfacesSize()

int Cell::getCellInterfacesSize ( ) const

◆ getChildVector()

std::vector< Cell * > * Cell::getChildVector ( )

Return pointer to child vector

◆ getCons()

Flux * Cell::getCons ( ) const

◆ getConsTransport()

Transport * Cell::getConsTransport ( const int &  numTransport) const

◆ getDataToReceiveAndRefine()

void Cell::getDataToReceiveAndRefine ( std::vector< double > &  dataToReceive,
std::vector< int > &  dataSplitToReceive,
const int &  lvl,
Eos **  eos,
int &  counter,
int &  counterSplit,
const int &  nbCellsY,
const int &  nbCellsZ,
const std::vector< AddPhys * > &  addPhys 
)

◆ getDensityGradient()

double Cell::getDensityGradient ( )

◆ getElement()

Element * Cell::getElement ( ) const

◆ getGradMixture()

virtual GradMixture * Cell::getGradMixture ( ) const
inlinevirtual

Reimplemented in CellO2, and CellO2NS.

◆ getGradPhase()

virtual GradPhase * Cell::getGradPhase ( const int &  ) const
inlinevirtual

Reimplemented in CellO2, and CellO2NS.

◆ getGradTk()

const Coord & Cell::getGradTk ( int &  numPhase,
int &  numAddPhys 
) const

Allow to recover an additional physical quantity.

◆ getGradTransport()

virtual GradTransport * Cell::getGradTransport ( const int &  ) const
inlinevirtual

Reimplemented in CellO2, and CellO2NS.

◆ getLvl()

const int & Cell::getLvl ( ) const
inline

◆ getMixture()

Mixture * Cell::getMixture ( Prim  = vecPhases) const
virtual

Reimplemented in CellO2.

◆ getModel()

Model * Cell::getModel ( )

◆ getNumberCellsChildren()

int Cell::getNumberCellsChildren ( )

Return the number of children cells

◆ getNumberPhases()

const int & Cell::getNumberPhases ( ) const

◆ getNumberTransports()

const int & Cell::getNumberTransports ( ) const

◆ getPhase()

Phase * Cell::getPhase ( const int &  phaseNumber,
Prim  = vecPhases 
) const
virtual

Reimplemented in CellO2.

◆ getPhases()

Phase ** Cell::getPhases ( Prim  = vecPhases) const
virtual

Reimplemented in CellO2.

◆ getPosition()

const Coord & Cell::getPosition ( ) const
inline

◆ getPsat()

double Cell::getPsat ( )

Compute saturation pressure for a liquid/vapor fluid mixture (first phase is considered predominant -> not generalized be careful)

◆ getQPA()

const QuantitiesAddPhys * Cell::getQPA ( int &  numQPA) const
inline

◆ getRankOfNeighborCPU()

virtual int Cell::getRankOfNeighborCPU ( ) const
inlinevirtual

Does nothing for non-ghost O2 cells

Reimplemented in CellGhost, CellO2GhostCartesian, and CellO2GhostNS.

◆ getSize()

const Coord & Cell::getSize ( ) const
inline

◆ getSizeX()

const double & Cell::getSizeX ( ) const
inline

◆ getSizeY()

const double & Cell::getSizeY ( ) const
inline

◆ getSizeZ()

const double & Cell::getSizeZ ( ) const
inline

◆ getSlopes()

virtual Phase * Cell::getSlopes ( ) const
inlinevirtual

Does nothing for first order cells

◆ getSlopesTransport()

virtual Transport * Cell::getSlopesTransport ( ) const
inlinevirtual

Does nothing for first order cells

◆ getSplit()

const bool & Cell::getSplit ( ) const
inline

Get the cell AMR level in the AMR tree

◆ getTransport()

Transport & Cell::getTransport ( const int &  numTransport,
Prim  = vecPhases 
) const
virtual

Reimplemented in CellO2.

◆ getTransports()

Transport * Cell::getTransports ( Prim  = vecPhases) const
virtual

Reimplemented in CellO2.

◆ getVecQuantitiesAddPhys()

std::vector< QuantitiesAddPhys * > & Cell::getVecQuantitiesAddPhys ( )

◆ getVelocity() [1/2]

Coord & Cell::getVelocity ( )

◆ getVelocity() [2/2]

const Coord & Cell::getVelocity ( ) const

◆ getWall()

bool Cell::getWall ( ) const
inline

◆ getXi()

const double & Cell::getXi ( ) const
inline

Return true if the cells is plit, false otherwise

◆ hasNeighboringGhostCellOfCPUneighbour()

bool Cell::hasNeighboringGhostCellOfCPUneighbour ( const int &  neighbour) const

Return a bool that is true if the cell has a neighboring ghost cell corresponding to CPU "neighbour"

◆ isCellGhost()

virtual bool Cell::isCellGhost ( ) const
inlinevirtual

Does nothing for first order cells

Reimplemented in CellGhost, and CellO2GhostCartesian.

◆ limitGradientsO2()

virtual void Cell::limitGradientsO2 ( Limiter )
inlinevirtual

Reimplemented in CellO2, and CellO2NS.

◆ localProjection()

void Cell::localProjection ( const Coord normal,
const Coord tangent,
const Coord binormal,
Prim  = vecPhases 
)

◆ lookForPmax()

void Cell::lookForPmax ( double *  pMax,
double *  pMaxWall 
)

◆ lvlNeighborTooHigh()

bool Cell::lvlNeighborTooHigh ( )

Look for AMR level of neighboring cells if too high to unrefine

◆ lvlNeighborTooLow()

bool Cell::lvlNeighborTooLow ( )

Look for AMR level of neighboring cells if too low to refine

◆ numberOfNeighboringGhostCellsOfCPUneighbour()

int Cell::numberOfNeighboringGhostCellsOfCPUneighbour ( const int &  neighbour) const

Return the number of neighboring ghost cells corresponding to CPU "neighbour" this cell has

◆ popBackSlope()

virtual void Cell::popBackSlope ( )
inlinevirtual

Does nothing for non-ghost O2 cells

Reimplemented in CellO2GhostCartesian.

◆ predictionOrdre2()

virtual void Cell::predictionOrdre2 ( const double &  ,
Symmetry  
)
inlinevirtual

Does nothing for first order cells

Reimplemented in CellO2.

◆ prepareAddPhys()

void Cell::prepareAddPhys ( )

◆ printGnuplotAMR()

bool Cell::printGnuplotAMR ( std::ofstream &  fileStream,
const int &  dim,
GeometricObject objet = 0,
bool  recordPsat = false 
)

◆ printInfo()

void Cell::printInfo ( ) const

◆ printPhasesMixture()

void Cell::printPhasesMixture ( std::ofstream &  fileStream) const

◆ pushBackSlope()

virtual void Cell::pushBackSlope ( )
inlinevirtual

Reimplemented in CellO2GhostCartesian.

◆ refineCellAndCellInterfaces()

void Cell::refineCellAndCellInterfaces ( const int &  nbCellsY,
const int &  nbCellsZ,
const std::vector< AddPhys * > &  addPhys,
const bool &  refineExternalCellInterfaces 
)

Refine parent cell by creation of children cells

◆ refineCellAndCellInterfacesGhost()

void Cell::refineCellAndCellInterfacesGhost ( const int &  nbCellsY,
const int &  nbCellsZ,
const std::vector< AddPhys * > &  addPhys 
)

Refinement of parent ghost cell by creation of children ghost cells

◆ reinitializeColorFunction()

void Cell::reinitializeColorFunction ( const int &  numTransport,
const int &  numPhase 
)

Re-initialize the color function (transport) with alpha.

◆ reverseProjection()

void Cell::reverseProjection ( const Coord normal,
const Coord tangent,
const Coord binormal 
)
virtual

◆ saveCons()

virtual void Cell::saveCons ( )
inlinevirtual

Does nothing for first order cells

Reimplemented in CellO2.

◆ selectScalar()

double Cell::selectScalar ( Variable  nameVariable,
int  num = 0 
) const

Select a specific scalar variable.

Parameters
nameVariablesName of the variable to select
numPhasesPhases number's

◆ selectVector()

Coord Cell::selectVector ( Variable  nameVector,
int  num = 0,
int  subscript = -1 
) const

◆ setCellInterface()

void Cell::setCellInterface ( const int &  b,
CellInterface cellInterface 
)

◆ setCons()

void Cell::setCons ( Flux cons)

◆ setConsTransport()

void Cell::setConsTransport ( double  value,
const int &  numTransport 
)

◆ setElement()

void Cell::setElement ( Element element,
const int &  numCell 
)

◆ setGradTk()

void Cell::setGradTk ( int &  numPhase,
int &  numAddPhys,
double *  buffer,
int &  counter 
)

◆ setRankOfNeighborCPU()

virtual void Cell::setRankOfNeighborCPU ( int  )
inlinevirtual

Reimplemented in CellGhost, CellO2GhostCartesian, and CellO2GhostNS.

◆ setToZeroCons()

void Cell::setToZeroCons ( )

◆ setToZeroConsGlobal()

void Cell::setToZeroConsGlobal ( )

◆ setToZeroConsXi()

void Cell::setToZeroConsXi ( )

set m_consXi to zero

◆ setToZeroXi()

void Cell::setToZeroXi ( )

set m_xi to zero

◆ setTransport()

void Cell::setTransport ( double  value,
int &  numTransport,
Prim  = vecPhases 
)
virtual

Reimplemented in CellO2.

◆ setVector()

void Cell::setVector ( Variable  nameVector,
const Coord value,
int  num = 0,
int  subscript = -1 
)

◆ setWall()

void Cell::setWall ( bool  wall)

◆ setXi()

void Cell::setXi ( double  value)

Return Xi cell value Set the Xi cell value

◆ sourceTermIntegration()

void Cell::sourceTermIntegration ( const double &  )
inline

◆ timeEvolution()

void Cell::timeEvolution ( const double &  dt,
Symmetry symmetry 
)

◆ timeEvolutionAddPhys()

void Cell::timeEvolutionAddPhys ( const double &  dt)

◆ timeEvolutionXi()

void Cell::timeEvolutionXi ( )

time evolution of Xi for smoothing

◆ traverseObjet()

bool Cell::traverseObjet ( const GeometricObject objet) const

◆ unrefineCellAndCellInterfaces()

void Cell::unrefineCellAndCellInterfaces ( )

Unrefine parent cell by destruction of children cells

◆ unrefineCellAndCellInterfacesGhost()

void Cell::unrefineCellAndCellInterfacesGhost ( )

Unrefinement of parent ghost cell by destruction of children ghost cells

◆ updateNbCellsTotalAMR()

void Cell::updateNbCellsTotalAMR ( int &  nbCellsTotalAMR)

◆ updatePointersInternalCellInterfaces()

void Cell::updatePointersInternalCellInterfaces ( )

Member Data Documentation

◆ m_cellInterfaces

std::vector<CellInterface*> Cell::m_cellInterfaces
protected

Vector of cell-interface pointers

◆ m_childrenCells

std::vector<Cell*> Cell::m_childrenCells
protected

Vector of children cells pointers

◆ m_childrenInternalCellInterfaces

std::vector<CellInterface*> Cell::m_childrenInternalCellInterfaces
protected

Vector of Internal children cell-interface pointers of the cell

◆ m_cons

Flux* Cell::m_cons
protected

Conservative variables

◆ m_consTransports

Transport* Cell::m_consTransports
protected

Array of fluxes for advected variables

◆ m_consXi

double Cell::m_consXi
protected

Buffer variable for Xi fluxes

◆ m_element

Element* Cell::m_element
protected

Pointer to corresponding geometrical mesh element

◆ m_lvl

int Cell::m_lvl
protected

Cell AMR level in the AMR tree

◆ m_mixture

Mixture* Cell::m_mixture
protected

◆ m_split

bool Cell::m_split
protected

Indicator for splitted cell (Do I possess children ?)

◆ m_vecPhases

Phase** Cell::m_vecPhases
protected

Array of phases

◆ m_vecQuantitiesAddPhys

std::vector<QuantitiesAddPhys*> Cell::m_vecQuantitiesAddPhys
protected

Vector of pointers to the Quantities of Additional Physics of the cell

◆ m_vecTransports

Transport* Cell::m_vecTransports
protected

Array of passive variables advected in the flow

◆ m_wall

bool Cell::m_wall
protected

Bool indicating if the cell is a solid boundary (for immersed boundaries)

◆ m_xi

double Cell::m_xi
protected

Criteria for refine/unrefine cell


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