CellGPU  0.8.0
GPU-accelerated simulations of cells
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
voronoiModelBase Class Reference

Perform and test triangulations in an MD setting, using kernels in voronoiModelBase Kernels. More...

#include <voronoiModelBase.h>

Inheritance diagram for voronoiModelBase:
Inheritance graph
[legend]

Public Member Functions

 voronoiModelBase ()
 The constructor! More...
 
void initializeVoronoiModelBase (int n)
 A default initialization scheme. More...
 
void setCPU (bool global=true)
 Enforce CPU-only operation. More...
 
void writeTriangulation (ofstream &outfile)
 write triangulation to text file
 
void readTriangulation (ifstream &infile)
 read positions from text file...for debugging
 
virtual void enforceTopology ()
 update/enforce the topology More...
 
void setExclusions (vector< int > &exes)
 Declare which particles are to be excluded (exes[i]!=0) More...
 
virtual int getNumberOfDegreesOfFreedom ()
 In voronoi models the number of degrees of freedom is the number of cells.
 
virtual void moveDegreesOfFreedom (GPUArray< Dscalar2 > &displacements, Dscalar scale=1.)
 moveDegrees of Freedom calls either the move points or move points CPU routines More...
 
virtual void getForces (GPUArray< Dscalar2 > &forces)
 return the forces
 
virtual GPUArray< Dscalar2 > & returnForces ()
 return a reference to the GPUArray of the current forces
 
virtual void computeGeometryCPU ()
 Compute cell geometry on the CPU. More...
 
virtual void computeGeometryGPU ()
 call gpu_compute_geometry kernel caller More...
 
virtual void cellDivision (const vector< int > &parameters, const vector< Dscalar > &dParams)
 allow for cell division, according to a vector of model-dependent parameters More...
 
virtual void cellDeath (int cellIndex)
 Kill the indexed cell by simply removing it from the simulation. More...
 
void movePoints (GPUArray< Dscalar2 > &displacements, Dscalar scale)
 move particles on the GPU More...
 
void movePointsCPU (GPUArray< Dscalar2 > &displacements, Dscalar scale)
 move particles on the CPU More...
 
void resetDelLocPoints ()
 Transfer particle data from the GPU to the CPU for use by delLoc. More...
 
void updateCellList ()
 Update the cell list structure after particles have moved. More...
 
void updateNeighIdxs ()
 update the NieghIdxs data structure More...
 
void initializeSimple2DActiveCell (int n)
 initialize class' data structures and set default values More...
 
void setv0Dr (Dscalar v0new, Dscalar drnew)
 Set uniform motility. More...
 
void setCellMotility (vector< Dscalar > &v0s, vector< Dscalar > &drs)
 Set non-uniform cell motilites. More...
 
void setCellDirectorsRandomly ()
 Set random cell directors (for active cell models) More...
 
Dscalar vicsekOrderParameter (Dscalar2 &vParallel, Dscalar2 &vPerpendicular)
 measure the viscek order parameter N^-1 {v_i}{|v_i}
 
Dscalar vicsekOrderParameterDirector (Dscalar2 &vParallel, Dscalar2 &vPerpendicular)
 measure the viscek order parameter N^-1 {v_i}{|v_i} from the director only
 
void initializeSimple2DCell (int n)
 initialize class' data structures and set default values More...
 
virtual void setGPU ()
 Enforce GPU-only operation. This is the default mode, so this method need not be called most of the time.
 
virtual void setCPU ()
 Enforce CPU-only operation. Derived classes might have to do more work when the CPU mode is invoked.
 
virtual void computeForces ()
 do everything necessary to compute forces in the current model
 
virtual void computeGeometry ()
 call either the computeGeometryCPU or GPU routines for the current model More...
 
virtual Dscalar computeEnergy ()
 do everything necessary to compute the energy for the current model
 
Dscalar computeKineticEnergy ()
 Call masses and velocities to get the total kinetic energy. More...
 
Dscalar4 computeKineticPressure ()
 Call masses and velocities to get the average kinetic contribution to the pressure tensor. More...
 
void setCellPreferencesUniform (Dscalar A0, Dscalar P0)
 Set uniform cell area and perimeter preferences. More...
 
void setCellPreferences (vector< Dscalar2 > &AreaPeriPreferences)
 Set cell area and perimeter preferences according to input vector. More...
 
void setCellPositionsRandomly ()
 Set random cell positions, and set the periodic box to a square with average cell area=1. More...
 
void setCellPositions (vector< Dscalar2 > newCellPositions)
 Set cell positions according to a user-specified vector. More...
 
void setVertexPositions (vector< Dscalar2 > newVertexPositions)
 Set vertex positions according to a user-specified vector. More...
 
Dscalar setCellVelocitiesMaxwellBoltzmann (Dscalar T)
 Set velocities via a temperature. The return value is the total kinetic energy. More...
 
Dscalar setVertexVelocitiesMaxwellBoltzmann (Dscalar T)
 Set velocities via a temperature for the vertex degrees of freedom. More...
 
void setModuliUniform (Dscalar newKA, Dscalar newKP)
 set uniform moduli for all cells More...
 
void setCellTypeUniform (int i)
 Set all cells to the same "type". More...
 
void setCellType (vector< int > &types)
 Set cells to different "type". More...
 
void setVertexTopologyFromCells (vector< vector< int > > cellVertexIndices)
 An uncomfortable function to allow the user to set vertex topology "by hand". More...
 
virtual gpuboxreturnBox ()
 return the gpubox
 
void setBox (BoxPtr _box)
 This can be used, but should not normally be. This re-assigns the pointer.
 
virtual vector< int > & returnItt ()
 return the base "itt" re-indexing vector
 
virtual GPUArray< Dscalar2 > & returnModuli ()
 Return a reference to moduli.
 
virtual GPUArray< Dscalar2 > & returnAreaPeri ()
 Return a reference to AreaPeri array.
 
virtual GPUArray< Dscalar2 > & returnAreaPeriPreferences ()
 Return a reference to AreaPeriPreferences.
 
virtual GPUArray< Dscalar2 > & returnVelocities ()
 Return a reference to velocities on cells. VertexModelBase will instead return vertexVelocities.
 
virtual GPUArray< Dscalar2 > & returnPositions ()
 Return a reference to Positions on cells.
 
virtual GPUArray< Dscalar > & returnMasses ()
 Return a reference to Masses on cells.
 
virtual GPUArray< Dscalar > & returnOtherData ()
 Return other data just returns the masses; in this class it's not needed.
 
void setDeltaT (Dscalar dt)
 Set the simulation time stepsize.
 
void getCellNeighs (int idx, int &nNeighs, vector< int > &neighs)
 
Dscalar getMaxForce ()
 Get the maximum force on a cell.
 
void reportMeanCellForce (bool verbose)
 Report the current average force on each cell. More...
 
void reportMeanVertexForce (bool verbose=false)
 Report the current average force per vertex...should be close to zero.
 
void reportAP (bool verbose=false)
 report the current total area, and optionally the area and perimeter for each cell
 
Dscalar reportq ()
 Report the average value of p/sqrt(A) for the cells in the system. More...
 
Dscalar reportVarq ()
 Report the variance of p/sqrt(A) for the cells in the system. More...
 
Dscalar2 reportVarAP ()
 Report the variance of A and P for the cells in the system. More...
 
Dscalar reportMeanP ()
 Report the mean value of the perimeter. More...
 
virtual void getDynMatEntries (vector< int2 > &rcs, vector< Dscalar > &vals, Dscalar unstress=1.0, Dscalar stress=1.0)
 Do whatever is necessary to get lists of dynamical matrix elements.
 
virtual void setTime (Dscalar time)
 set the time
 

Public Attributes

DelaunayLoc delLoc
 The class' local Delaunay tester/updater.
 
Dscalar repPerFrame
 Collect statistics of how many triangulation repairs are done per frame, etc.
 
int skippedFrames
 How often were all circumcenters empty (so that no data transfers and no repairs were necessary)?
 
int GlobalFixes
 How often were global re-triangulations performed?
 
GPUArray< Dscalar2 > external_forces
 "exclusions" zero out the force on a cell...the external force needed to do this is stored in external_forces
 
GPUArray< int > exclusions
 An array containing the indices of excluded particles.
 
int localTopologyUpdates
 The number of topology updates performed at the individual particle level.
 
GPUArray< Dscalar > cellDirectors
 An array of angles (relative to the x-axis) that the cell directors point.
 
GPUArray< Dscalar > cellDirectorForces
 An array of forces acting on the cell directors.
 
Dscalar v0
 velocity of cells in mono-motile systems
 
Dscalar Dr
 rotational diffusion of cell directors in mono-motile systems
 
GPUArray< Dscalar2 > Motility
 The motility parameters (v0 and Dr) for each cell.
 
int Ncells
 Number of cells in the simulation.
 
int Nvertices
 Number of vertices.
 
GPUArray< Dscalar2 > cellPositions
 Cell positions... not used for computation, but can track, e.g., MSD of cell centers.
 
GPUArray< Dscalar2 > vertexPositions
 Position of the vertices.
 
GPUArray< Dscalar2 > cellVelocities
 The velocity vector of cells (only relevant if the equations of motion use it)
 
GPUArray< Dscalar > cellMasses
 The masses of the cells.
 
GPUArray< Dscalar2 > vertexVelocities
 The velocity vector of vertices (only relevant if the equations of motion use it)
 
GPUArray< Dscalar > vertexMasses
 The masses of the vertices.
 
GPUArray< int > vertexNeighbors
 VERTEX neighbors of every vertex. More...
 
GPUArray< int > vertexCellNeighbors
 Cell neighbors of every vertex. More...
 
Index2D n_idx
 A 2dIndexer for computing where in the GPUArray to look for a given cell's vertices.
 
GPUArray< int > cellNeighborNum
 The number of CELL neighbors of each cell. For simple models this is the same as cellVertexNum, but does not have to be.
 
GPUArray< int > cellNeighbors
 CELL neighbors of every cell.
 
GPUArray< int > cellVertexNum
 The number of vertices defining each cell. More...
 
GPUArray< Dscalar2 > vertexForces
 an array containing net force on each vertex
 
GPUArray< Dscalar2 > cellForces
 an array containing net force on each cell
 
GPUArray< int > cellType
 An array of integers labeling cell type...an easy way of determining if cells are different. More...
 
Index2D cellTypeIndexer
 A indexer for turning a pair of cells into a 1-D index.
 
Dscalar Energy
 The current potential energy of the system; only updated when an explicit energy calculation is called (i.e. not by default each timestep)
 
Dscalar KineticEnergy
 The current kinetic energy of the system; only updated when an explicit calculation is called.
 
vector< int > tagToIdx
 To write consistent files...the cell that started the simulation as index i has current index tagToIdx[i]. More...
 
vector< int > tagToIdxVertex
 To write consistent files...the vertex that started the simulation as index i has current index tagToIdx[i].
 
BoxPtr Box
 the box defining the periodic domain
 
int Timestep
 Count the number of times "performTimeStep" has been called.
 
Dscalar deltaT
 The time stepsize of the simulation.
 
bool forcesUpToDate
 Are the forces (and hence, the geometry) up-to-date?
 
Dscalar currentTime
 a time variable for keeping track of the simulation variable (for databases)
 

Protected Member Functions

void spatialSorting ()
 sort points along a Hilbert curve for data locality More...
 
void fullTriangulation (bool verbose=false)
 construct the global periodic triangulation point-by-point using non-CGAL methods More...
 
void globalTriangulationCGAL (bool verbose=false)
 Globally construct the triangulation via CGAL. More...
 
void getCircumcenterIndices (bool secondtime=false, bool verbose=false)
 build the auxiliary data structure containing the indices of the particle circumcenters from the neighbor list More...
 
void testTriangulation ()
 Test the current neighbor list to see if it is still a valid triangulation. GPU function. More...
 
void testTriangulationCPU ()
 Test the validity of the triangulation on the CPU. More...
 
void repairTriangulation (vector< int > &fixlist)
 repair any problems with the triangulation on the CPU More...
 
void testAndRepairTriangulation (bool verb=false)
 A workhorse function that calls the appropriate topology testing and repairing routines. More...
 
void allDelSets ()
 call getDelSets for all particles More...
 
bool getDelSets (int i)
 
void resetLists ()
 resize all neighMax-related arrays More...
 
void resizeAndReset ()
 do resize and resetting operations common to cellDivision and cellDeath More...
 
Matrix2x2 dHdri (Dscalar2 ri, Dscalar2 rj, Dscalar2 rk)
 The derivative of a voronoi vertex position with respect to change in the first cells position. More...
 
Dscalar2 dAidrj (int i, int j)
 Derivative of the area of cell i with respect to the position of cell j. More...
 
Dscalar2 dPidrj (int i, int j)
 Derivative of the perimeter of cell i with respect to the position of cell j. More...
 
Matrix2x2 d2Areadvdr (Matrix2x2 &dvpdr, Matrix2x2 &dvmdr)
 Second derivative of area w/r/t voronoi and cell position. More...
 
Matrix2x2 d2Peridvdr (Matrix2x2 &dvdr, Matrix2x2 &dvmdr, Matrix2x2 &dvpdr, Dscalar2 vm, Dscalar2 v, Dscalar2 vp)
 Second derivative of perimeter w/r/t voronoi and cell position. More...
 
vector< Dscalar > d2Hdridrj (Dscalar2 rj, Dscalar2 rk, int jj)
 second derivatives of voronoi vertex with respect to cell positions More...
 
void spatiallySortVerticesAndCellActivity ()
 call the Simple2DCell spatial vertex sorter, and re-index arrays of cell activity More...
 
void spatiallySortCellsAndCellActivity ()
 call the Simple2DCell spatial cell sorter, and re-index arrays of cell activity More...
 
void initializeCellSorting ()
 set the size of the cell-sorting structures, initialize lists simply More...
 
void initializeVertexSorting ()
 set the size of the vertex-sorting structures, initialize lists simply More...
 
void reIndexCellArray (GPUArray< int > &array)
 Re-index cell arrays after a spatial sorting has occured. More...
 
void reIndexCellArray (GPUArray< Dscalar > &array)
 why use templates when you can type more? More...
 
void reIndexCellArray (GPUArray< Dscalar2 > &array)
 why use templates when you can type more? More...
 
void reIndexVertexArray (GPUArray< int > &array)
 Re-index vertex after a spatial sorting has occured.
 
void reIndexVertexArray (GPUArray< Dscalar > &array)
 why use templates when you can type more?
 
void reIndexVertexArray (GPUArray< Dscalar2 > &array)
 why use templates when you can type more? More...
 
void spatiallySortCells ()
 Perform a spatial sorting of the cells to try to maintain data locality. More...
 
void spatiallySortVertices ()
 Perform a spatial sorting of the vertices to try to maintain data locality. More...
 

Protected Attributes

cellListGPU celllist
 The associated cell list structure.
 
Dscalar cellsize
 The size of the cell list's underlying grid.
 
int neighMax
 An upper bound for the maximum number of neighbors that any cell has.
 
GPUArray< int2 > NeighIdxs
 
int NeighIdxNum
 A utility integer to help with NeighIdxs.
 
GPUArray< int3 > circumcenters
 A data structure that holds the indices of particles forming the circumcircles of the Delaunay Triangulation.
 
int NumCircumCenters
 The number of circumcircles...for a periodic system, this should never change. This provides one check that local updates to the triangulation are globally consistent.
 
GPUArray< int > anyCircumcenterTestFailed
 A flag that can be accessed by child classes... serves as notification that any change in the network topology has occured.
 
int completeRetriangulationPerformed
 A flag that notifies that a global re-triangulation has been performed.
 
bool neighMaxChange
 A flag that notifies that the maximum number of neighbors may have changed, necessitating resizing of some data arrays.
 
GPUArray< int > repair
 A a vector of zeros (everything is fine) and ones (that index needs to be repaired)
 
vector< int > NeedsFixing
 A smaller vector that, after testing the triangulation, contains the particle indices that need their local topology to be updated.
 
bool globalOnly
 When true, the CPU branch will execute global retriangulations through CGAL on every time step. More...
 
int timestep
 Count the number of times that testAndRepair has been called, separately from the derived class' time.
 
bool particleExclusions
 A flag that notifies the existence of any particle exclusions (for which the net force is set to zero by fictitious external forces)
 
GPUArray< int2 > delSets
 delSet.data[n_idx(nn,i)] are the previous and next consecutive delaunay neighbors More...
 
GPUArray< int > delOther
 delOther.data[n_idx(nn,i)] contains the index of the "other" delaunay neighbor. More...
 
GPUArray< Dscalar2 > forceSets
 In GPU mode, interactions are computed "per voronoi vertex"...forceSets are summed up to get total force on a particle.
 
bool GPUcompute
 Compute aspects of the model on the GPU.
 
bool Reproducible
 A flag that determines whether the GPU RNG is the same every time.
 
noiseSource noise
 A source of noise for random cell initialization.
 
Dscalar KA
 the area modulus
 
Dscalar KP
 The perimeter modulus.
 
GPUArray< Dscalar2 > Moduli
 The area and perimeter moduli of each cell. CURRENTLY NOT SUPPORTED, BUT EASY TO IMPLEMENT.
 
GPUArray< Dscalar2 > AreaPeri
 The current area and perimeter of each cell.
 
GPUArray< Dscalar2 > AreaPeriPreferences
 The area and perimeter preferences of each cell.
 
GPUArray< int > cellVertices
 A structure that indexes the vertices defining each cell. More...
 
int vertexMax
 An upper bound for the maximum number of neighbors that any cell has.
 
GPUArray< Dscalar2 > voroCur
 3*Nvertices length array of the position of vertices around cells More...
 
GPUArray< Dscalar4 > voroLastNext
 
vector< int > itt
 A map between cell index and the spatially sorted version. More...
 
vector< int > tti
 A temporary structure that inverts itt.
 
vector< int > idxToTag
 A temporary structure that inverse tagToIdx.
 
vector< int > ittVertex
 A map between vertex index and the spatially sorted version.
 
vector< int > ttiVertex
 A temporary structure that inverts itt.
 
vector< int > idxToTagVertex
 A temporary structure that inverse tagToIdx.
 
GPUArray< Dscalar2 > displacements
 An array of displacements used only for the equations of motion.
 

Detailed Description

Perform and test triangulations in an MD setting, using kernels in voronoiModelBase Kernels.

voronoiModelBase is a core engine class, capable of taking a set of points in a periodic domain, performing Delaunay triangulations on them, testing whether those triangulations are valid on either the CPU or GPU, and locally repair invalid triangulations on the CPU.

Voronoi models have their topology taken care of by the underlying triangulation, and so child classes just need to implement an energy functions (and corresponding force law)

Constructor & Destructor Documentation

◆ voronoiModelBase()

voronoiModelBase::voronoiModelBase ( )

The constructor!

A simple constructor that sets many of the class member variables to zero

Member Function Documentation

◆ initializeVoronoiModelBase()

void voronoiModelBase::initializeVoronoiModelBase ( int  n)

A default initialization scheme.

a function that takes care of the initialization of the class.

Parameters
nthe number of cells to initialize

◆ setCPU()

void voronoiModelBase::setCPU ( bool  global = true)
inlinevirtual

Enforce CPU-only operation.

Parameters
globaldefaults to true. When global is set to true, the CPU branch will try the local repair scheme. This is generally slower, but if working in a regime where things change very infrequently, it may be faster.

Implements Simple2DCell.

References enforceTopology(), globalOnly, Simple2DCell::GPUcompute, readTriangulation(), setExclusions(), and writeTriangulation().

◆ enforceTopology()

void voronoiModelBase::enforceTopology ( )
virtual

update/enforce the topology

goes through the process of testing and repairing the topology on either the CPU or GPU

Postcondition
and topological changes needed by cell motion are detected and repaired

Reimplemented from Simple2DCell.

References allDelSets(), anyCircumcenterTestFailed, completeRetriangulationPerformed, ArrayHandle< T >::data, delOther, delSets, access_location::device, getDelSets(), globalTriangulationCGAL(), Simple2DCell::GPUcompute, access_location::host, NeedsFixing, NeighIdxs, neighMaxChange, access_mode::read, resetLists(), and testAndRepairTriangulation().

Referenced by setCPU().

◆ setExclusions()

void voronoiModelBase::setExclusions ( vector< int > &  exes)

Declare which particles are to be excluded (exes[i]!=0)

Parameters
exesa list of per-particle indications of whether a particle should be excluded (exes[i] !=0) or not/

References ArrayHandle< T >::data, exclusions, external_forces, access_location::host, Simple2DActiveCell::Motility, Simple2DCell::Ncells, access_mode::overwrite, particleExclusions, access_mode::readwrite, and GPUArray< T >::resize().

Referenced by setCPU().

◆ moveDegreesOfFreedom()

void voronoiModelBase::moveDegreesOfFreedom ( GPUArray< Dscalar2 > &  displacements,
Dscalar  scale = 1. 
)
virtual

moveDegrees of Freedom calls either the move points or move points CPU routines

Displace cells on either the GPU or CPU, according to the flag

Parameters
displacementsa vector of Dscalar2 specifying how much to move every cell
scalea scalar that multiples the value of every index of displacements before things are moved
Postcondition
the cells are displaced according to the input vector, and then put back in the main unit cell.

Reimplemented from Simple2DCell.

Referenced by getNumberOfDegreesOfFreedom().

◆ computeGeometryCPU()

void voronoiModelBase::computeGeometryCPU ( )
virtual

Compute cell geometry on the CPU.

Precondition
Topology is up-to-date on the CPU
Postcondition
geometry and voronoi neighbor locations are computed for the current configuration

Reimplemented from Simple2DCell.

References Simple2DCell::AreaPeri, Simple2DCell::cellNeighborNum, Simple2DCell::cellNeighbors, Simple2DCell::cellPositions, ArrayHandle< T >::data, access_location::host, Simple2DCell::n_idx, Simple2DCell::Ncells, access_mode::overwrite, access_mode::read, access_mode::readwrite, Simple2DCell::voroCur, and Simple2DCell::voroLastNext.

Referenced by returnForces().

◆ computeGeometryGPU()

void voronoiModelBase::computeGeometryGPU ( )
virtual

call gpu_compute_geometry kernel caller

Precondition
The topology of the Delaunay triangulation is up-to-date on the GPU
Postcondition
calculate all cell areas, perimenters, and voronoi neighbors

Reimplemented from Simple2DCell.

References Simple2DCell::AreaPeri, Simple2DCell::Box, Simple2DCell::cellNeighborNum, Simple2DCell::cellNeighbors, Simple2DCell::cellPositions, ArrayHandle< T >::data, access_location::device, gpu_compute_voronoi_geometry(), Simple2DCell::n_idx, Simple2DCell::Ncells, access_mode::overwrite, access_mode::read, access_mode::readwrite, Simple2DCell::voroCur, and Simple2DCell::voroLastNext.

Referenced by returnForces().

◆ cellDivision()

void voronoiModelBase::cellDivision ( const vector< int > &  parameters,
const vector< Dscalar > &  dParams 
)
virtual

allow for cell division, according to a vector of model-dependent parameters

Trigger a cell division event, which involves some laborious re-indexing of various data structures. This version uses a heavy-handed approach, hitting the cell positions with a full, global retriangulation (rather than just updating the targeted cell and its neighbor with the local topology repair routines). If a simulation is needed where the cell division rate is very rapid, this should be improved. The idea of the division is that a targeted cell will divide normal to an axis specified by the angle, theta, passed to the function. The final state cell positions are placed along the axis at a distance away from the initial cell position set by a multiplicative factor (<1) of the in-routine determined maximum distance in the cell along that axis. parameters[0] = the index of the cell to undergo a division event dParams[0] = an angle, theta dParams[1] = a fraction of maximum separation of the new cell positions along the axis of the cell specified by theta

This function is meant to be called before the start of a new timestep. It should be immediately followed by a computeGeometry call.

Postcondition
the new cell is the final indexed entry of the various data structures (e.g., cellPositions[new number of cells - 1])

Reimplemented from Simple2DActiveCell.

Referenced by returnForces().

◆ cellDeath()

void voronoiModelBase::cellDeath ( int  cellIndex)
virtual

Kill the indexed cell by simply removing it from the simulation.

Trigger a cell death event. In the Voronoi model this simply removes the targeted cell and instantaneously computes the new tesselation. Very violent, if the cell isn't already small.

Reimplemented from Simple2DActiveCell.

References Simple2DActiveCell::cellDeath(), and resizeAndReset().

Referenced by returnForces().

◆ movePoints()

void voronoiModelBase::movePoints ( GPUArray< Dscalar2 > &  displacements,
Dscalar  scale 
)

move particles on the GPU

Displace cells on the GPU

Parameters
displacementsa vector of Dscalar2 specifying how much to move every cell
Postcondition
the cells are displaced according to the input vector, and then put back in the main unit cell.

Referenced by returnForces().

◆ movePointsCPU()

void voronoiModelBase::movePointsCPU ( GPUArray< Dscalar2 > &  displacements,
Dscalar  scale 
)

move particles on the CPU

Displace cells on the CPU

Parameters
displacementsa vector of Dscalar2 specifying how much to move every cell
Postcondition
the cells are displaced according to the input vector, and then put back in the main unit cell.

◆ resetDelLocPoints()

void voronoiModelBase::resetDelLocPoints ( )

Transfer particle data from the GPU to the CPU for use by delLoc.

The GPU moves the location of points in the GPU memory... this gets a local copy that can be used by the DelaunayLoc class

Postcondition
the DelaunayLoc class has the updated set of cell positions, and its cell list is initialized

References Simple2DCell::cellPositions, cellsize, delLoc, access_location::host, DelaunayLoc::initialize(), Simple2DCell::Ncells, access_mode::read, and DelaunayLoc::setPoints().

Referenced by fullTriangulation(), resizeAndReset(), and testTriangulationCPU().

◆ updateCellList()

void voronoiModelBase::updateCellList ( )

Update the cell list structure after particles have moved.

Postcondition
the cell list is updated according to the current cell positions

References celllist, Simple2DCell::cellPositions, cellListGPU::compute(), cellListGPU::computeGPU(), ArrayHandle< T >::data, Simple2DCell::GPUcompute, access_location::host, Simple2DCell::Ncells, access_mode::read, and cellListGPU::setParticles().

Referenced by testTriangulation().

◆ updateNeighIdxs()

void voronoiModelBase::updateNeighIdxs ( )

update the NieghIdxs data structure

Postcondition
the NeighIdx data structure is updated, which helps cut down on the number of inactive threads in the force set computation function

References Simple2DCell::cellNeighborNum, access_location::host, NeighIdxs, access_mode::overwrite, and access_mode::read.

Referenced by allDelSets(), and fullTriangulation().

◆ spatialSorting()

void voronoiModelBase::spatialSorting ( )
protectedvirtual

sort points along a Hilbert curve for data locality

When sortPeriod < 0, this routine does not get called

Postcondition
call Simple2DActiveCell's underlying Hilbert sort scheme, and re-index voronoiModelBase's extra arrays

Reimplemented from Simple2DModel.

References allDelSets(), exclusions, globalTriangulationCGAL(), Simple2DCell::reIndexCellArray(), resetLists(), and Simple2DActiveCell::spatiallySortCellsAndCellActivity().

◆ fullTriangulation()

void voronoiModelBase::fullTriangulation ( bool  verbose = false)
protected

◆ globalTriangulationCGAL()

void voronoiModelBase::globalTriangulationCGAL ( bool  verbose = false)
protected

Globally construct the triangulation via CGAL.

This function calls the DelaunayCGAL class to determine the Delaunay triangulation of the entire square periodic domain this method is, obviously, better than the version written by DMS, so should be the default option. In addition to performing a triangulation, the function also automatically calls updateNeighIdxs and getCircumcenterIndices/

References Simple2DCell::cellPositions, completeRetriangulationPerformed, GlobalFixes, access_location::host, Simple2DCell::Ncells, and access_mode::read.

Referenced by enforceTopology(), getCircumcenterIndices(), resizeAndReset(), spatialSorting(), testAndRepairTriangulation(), and testTriangulationCPU().

◆ getCircumcenterIndices()

void voronoiModelBase::getCircumcenterIndices ( bool  secondtime = false,
bool  verbose = false 
)
protected

build the auxiliary data structure containing the indices of the particle circumcenters from the neighbor list

Converts the neighbor list data structure into a list of the three particle indices defining all of the circumcenters in the triangulation. Keeping this version of the topology on the GPU allows for fast testing of what points need to be retriangulated.

References Simple2DCell::cellNeighborNum, Simple2DCell::cellNeighbors, circumcenters, ArrayHandle< T >::data, globalTriangulationCGAL(), access_location::host, Simple2DCell::n_idx, Simple2DCell::Ncells, NumCircumCenters, access_mode::overwrite, access_mode::read, timestep, and writeTriangulation().

Referenced by fullTriangulation().

◆ testTriangulation()

void voronoiModelBase::testTriangulation ( )
protected

Test the current neighbor list to see if it is still a valid triangulation. GPU function.

Call the GPU to test each circumcenter to see if it is still empty (i.e., how much of the triangulation from the last time step is still valid?). Note that because gpu_test_circumcenters always* copies at least a single integer back and forth (to answer the question "did any circumcircle come back non-empty?" for the cpu) this function is always an implicit cuda synchronization event. At least until non-default streams are added to the code.

References anyCircumcenterTestFailed, Simple2DCell::Box, cellListGPU::cell_indexer, cellListGPU::cell_list_indexer, cellListGPU::cell_sizes, celllist, Simple2DCell::cellPositions, circumcenters, ArrayHandle< T >::data, access_location::device, cellListGPU::getBoxsize(), cellListGPU::getXsize(), cellListGPU::getYsize(), gpu_test_circumcenters(), access_location::host, cellListGPU::idxs, Simple2DCell::Ncells, NumCircumCenters, access_mode::overwrite, access_mode::read, access_mode::readwrite, repair, and updateCellList().

Referenced by testAndRepairTriangulation().

◆ testTriangulationCPU()

void voronoiModelBase::testTriangulationCPU ( )
protected

Test the validity of the triangulation on the CPU.

perform the same check on the CPU... because of the cost of checking circumcircles and the relatively poor performance of the 1-ring calculation in DelaunayLoc, it is sometimes better to just re-triangulate the entire point set with CGAL. At the moment that is the default behavior of the cpu branch.

References anyCircumcenterTestFailed, Simple2DCell::cellNeighborNum, Simple2DCell::cellNeighbors, ArrayHandle< T >::data, delLoc, globalOnly, globalTriangulationCGAL(), access_location::host, localTopologyUpdates, Simple2DCell::n_idx, Simple2DCell::Ncells, neighMax, access_mode::readwrite, repair, resetDelLocPoints(), skippedFrames, and DelaunayLoc::testPointTriangulation().

Referenced by testAndRepairTriangulation().

◆ repairTriangulation()

void voronoiModelBase::repairTriangulation ( vector< int > &  fixlist)
protected

repair any problems with the triangulation on the CPU

Given a list of particle indices that need to be repaired, call CGAL to figure out their neighbors and then update the relevant data structures.

References repPerFrame.

Referenced by testAndRepairTriangulation().

◆ testAndRepairTriangulation()

void voronoiModelBase::testAndRepairTriangulation ( bool  verb = false)
protected

A workhorse function that calls the appropriate topology testing and repairing routines.

This function calls the relevant testing and repairing functions, and increments the "timestep" by one. Note that the call to testTriangulation will always synchronize the gpu (via a memcpy of the "anyCircumcenterTestFailed" variable)

References anyCircumcenterTestFailed, Simple2DCell::cellNeighborNum, Simple2DCell::cellNeighbors, completeRetriangulationPerformed, ArrayHandle< T >::data, globalTriangulationCGAL(), Simple2DCell::GPUcompute, access_location::host, Simple2DCell::n_idx, Simple2DCell::Ncells, NeedsFixing, access_mode::read, access_mode::readwrite, repair, repairTriangulation(), skippedFrames, testTriangulation(), testTriangulationCPU(), and timestep.

Referenced by enforceTopology().

◆ allDelSets()

void voronoiModelBase::allDelSets ( )
protected

call getDelSets for all particles

Calls updateNeighIdxs and then getDelSets(i) for all cells i

References getDelSets(), Simple2DCell::Ncells, and updateNeighIdxs().

Referenced by enforceTopology(), resizeAndReset(), and spatialSorting().

◆ getDelSets()

bool voronoiModelBase::getDelSets ( int  i)
protected

Maintain the delSets and delOther data structure for particle i If it returns false there was a problem and a global re-triangulation is triggered.

Parameters
ithe cell in question
Postcondition
the delSet and delOther data structure for cell i is updated. Recall that delSet.data[n_idx(nn,i)] is an int2; the x and y parts store the index of the previous and next Delaunay neighbor, ordered CCW. delOther contains the mutual neighbor of delSet.data[n_idx(nn,i)].y and delSet.data[n_idx(nn,i)].z that isn't cell i

References Simple2DCell::cellNeighborNum, Simple2DCell::cellNeighbors, ArrayHandle< T >::data, delOther, delSets, access_location::host, Simple2DCell::n_idx, access_mode::read, and access_mode::readwrite.

Referenced by allDelSets(), and enforceTopology().

◆ resetLists()

void voronoiModelBase::resetLists ( )
protected

resize all neighMax-related arrays

As the code is modified, all GPUArrays whose size depend on neighMax should be added to this function

Postcondition
voroCur,voroLastNext, delSets, delOther, and forceSets grow to size neighMax*Ncells

References delOther, delSets, forceSets, Simple2DCell::Ncells, neighMax, GPUArray< T >::resize(), Simple2DCell::voroCur, and Simple2DCell::voroLastNext.

Referenced by enforceTopology(), resizeAndReset(), and spatialSorting().

◆ resizeAndReset()

void voronoiModelBase::resizeAndReset ( )
protected

do resize and resetting operations common to cellDivision and cellDeath

A utility function for resizing data arrays... used by the cell division and cell death routines

References allDelSets(), Simple2DCell::cellForces, celllist, Simple2DCell::cellNeighborNum, circumcenters, Simple2DCell::displacements, exclusions, external_forces, globalTriangulationCGAL(), Simple2DCell::Ncells, NeighIdxs, repair, resetDelLocPoints(), resetLists(), GPUArray< T >::resize(), and cellListGPU::setNp().

Referenced by cellDeath().

◆ dHdri()

Matrix2x2 voronoiModelBase::dHdri ( Dscalar2  ri,
Dscalar2  rj,
Dscalar2  rk 
)
protected

The derivative of a voronoi vertex position with respect to change in the first cells position.

Parameters
riThe position of cell i
rjThe position of cell j
rkThe position of cell k Returns the derivative of the voronoi vertex shared by cells i, j , and k with respect to changing the position of cell i the (row, column) format specifies dH_{row}/dr_{i,column}

◆ dAidrj()

Dscalar2 voronoiModelBase::dAidrj ( int  i,
int  j 
)
protected

Derivative of the area of cell i with respect to the position of cell j.

Parameters
iThe index of cell i
jThe index of cell j
Precondition
Requires that computeGeometry is current Returns the derivative of the area of cell i w/r/t the position of cell j

◆ dPidrj()

Dscalar2 voronoiModelBase::dPidrj ( int  i,
int  j 
)
protected

Derivative of the perimeter of cell i with respect to the position of cell j.

Parameters
iThe index of cell i
jThe index of cell j Returns the derivative of the perimeter of cell i w/r/t the position of cell j

◆ d2Areadvdr()

Matrix2x2 voronoiModelBase::d2Areadvdr ( Matrix2x2 dvpdr,
Matrix2x2 dvmdr 
)
protected

Second derivative of area w/r/t voronoi and cell position.

Parameters
dvpdrderivative of v_{i+1} w/r/t a cell position
dvmdrderivative of v_{i-1} w/r/t a cell position

References Matrix2x2::x11.

◆ d2Peridvdr()

Matrix2x2 voronoiModelBase::d2Peridvdr ( Matrix2x2 dvdr,
Matrix2x2 dvmdr,
Matrix2x2 dvpdr,
Dscalar2  vm,
Dscalar2  v,
Dscalar2  vp 
)
protected

Second derivative of perimeter w/r/t voronoi and cell position.

Parameters
dvdrderivative of v_{i} w/r/t a cell position
dvpdrderivative of v_{i+1} w/r/t a cell position
dvmdrderivative of v_{i-1} w/r/t a cell position
vmposition of v_{i-1}
vposition of v_{i}
vpposition of v_{i+1}

◆ d2Hdridrj()

vector< Dscalar > voronoiModelBase::d2Hdridrj ( Dscalar2  rj,
Dscalar2  rk,
int  jj 
)
protected

second derivatives of voronoi vertex with respect to cell positions

Parameters
riThe position of cell i
rjThe position of cell j
rkThe position of cell k
jjthe index EITHER 1 or 2 of the second derivative Returns an 8-component vector containing the derivatives of the voronoi vertex formed by cells i, j, and k with respect to r_i and r_{jj}... jj should be either 1 (to give d^2H/(d r_i)^2 or 2 (to give d^2H/dridrj) The vector is laid out as (H_x/r_{i,x}r_{j,x}, H_y/r_{i,x}r_{j,x} H_x/r_{i,y}r_{j,x}, H_y/r_{i,y}r_{j,x H_x/r_{i,x}r_{j,y}, H_y/r_{i,x}r_{j,y H_x/r_{i,y}r_{j,y}, H_y/r_{i,y}r_{j,y} ) NOTE: This function does not check that ri, rj, and rk actually share a voronoi vertex in the triangulation NOTE: This function assumes that rj and rk are w/r/t the position of ri, so ri = (0.,0.)

◆ initializeSimple2DActiveCell()

void Simple2DActiveCell::initializeSimple2DActiveCell ( int  n)
inherited

initialize class' data structures and set default values

Initialize the data structures to the size specified by n, and set default values, and call Simple2DCell's initilization routine.

◆ setv0Dr()

void Simple2DActiveCell::setv0Dr ( Dscalar  v0new,
Dscalar  drnew 
)
virtualinherited

Set uniform motility.

Parameters
v0newthe new value of velocity for all cells
drnewthe new value of the rotational diffusion of cell directors for all cells

Implements Simple2DCell.

◆ setCellMotility()

void Simple2DActiveCell::setCellMotility ( vector< Dscalar > &  v0s,
vector< Dscalar > &  drs 
)
inherited

Set non-uniform cell motilites.

Parameters
v0sthe per-particle vector of what all velocities will be
drsthe per-particle vector of what all rotational diffusions will be

References Simple2DActiveCell::cellDirectors, Simple2DCell::cellVelocities, ArrayHandle< T >::data, access_location::host, Simple2DActiveCell::Motility, Simple2DCell::Ncells, access_mode::overwrite, access_mode::read, and GPUArray< T >::resize().

◆ setCellDirectorsRandomly()

void Simple2DActiveCell::setCellDirectorsRandomly ( )
inherited

◆ spatiallySortVerticesAndCellActivity()

void Simple2DActiveCell::spatiallySortVerticesAndCellActivity ( )
protectedinherited

call the Simple2DCell spatial vertex sorter, and re-index arrays of cell activity

Calls the spatial vertex sorting routine in Simple2DCell, and re-indexes the arrays for the cell RNGS, as well as the cell motility and cellDirector arrays

References Simple2DActiveCell::cellDirectors, Simple2DActiveCell::Motility, Simple2DCell::reIndexCellArray(), and Simple2DCell::spatiallySortVertices().

Referenced by vertexModelBase::spatialSorting().

◆ spatiallySortCellsAndCellActivity()

void Simple2DActiveCell::spatiallySortCellsAndCellActivity ( )
protectedinherited

call the Simple2DCell spatial cell sorter, and re-index arrays of cell activity

Calls the spatial vertex sorting routine in Simple2DCell, and re-indexes the arrays for the cell RNGS, as well as the cell motility and cellDirector arrays

References Simple2DActiveCell::cellDirectors, Simple2DActiveCell::Motility, Simple2DCell::reIndexCellArray(), and Simple2DCell::spatiallySortCells().

Referenced by spatialSorting().

◆ initializeSimple2DCell()

void Simple2DCell::initializeSimple2DCell ( int  n)
inherited

initialize class' data structures and set default values

Initialize the data structures to the size specified by n, and set default values.

◆ computeGeometry()

void Simple2DCell::computeGeometry ( )
virtualinherited

◆ computeKineticEnergy()

Dscalar Simple2DCell::computeKineticEnergy ( )
inherited

Call masses and velocities to get the total kinetic energy.

E = 0.5*m_i v_i^2

◆ computeKineticPressure()

Dscalar4 Simple2DCell::computeKineticPressure ( )
inherited

Call masses and velocities to get the average kinetic contribution to the pressure tensor.

P_ab = m_i v_{ib}v_{ia}

◆ setCellPreferencesUniform()

void Simple2DCell::setCellPreferencesUniform ( Dscalar  A0,
Dscalar  P0 
)
inherited

Set uniform cell area and perimeter preferences.

Generically believe that cells in 2D have a notion of a preferred area and perimeter

Referenced by Simple2DCell::enforceTopology().

◆ setCellPreferences()

void Simple2DCell::setCellPreferences ( vector< Dscalar2 > &  APPref)
inherited

Set cell area and perimeter preferences according to input vector.

Set the Area and Perimeter preferences to the input vector

References Simple2DCell::AreaPeriPreferences, ArrayHandle< T >::data, access_location::host, Simple2DCell::Ncells, access_mode::overwrite, and GPUArray< T >::resize().

◆ setCellPositionsRandomly()

void Simple2DCell::setCellPositionsRandomly ( )
inherited

Set random cell positions, and set the periodic box to a square with average cell area=1.

Resize the box so that every cell has, on average, area = 1, and place cells via either a simple, reproducible RNG or a non-reproducible RNG

References Simple2DCell::cellPositions, Simple2DCell::Ncells, and GPUArray< T >::resize().

◆ setCellPositions()

void Simple2DCell::setCellPositions ( vector< Dscalar2 >  newCellPositions)
inherited

Set cell positions according to a user-specified vector.

Does not update any other lists – it is the user's responsibility to maintain topology, etc, when using this function.

References Simple2DCell::cellPositions, ArrayHandle< T >::data, GPUArray< T >::getNumElements(), access_location::host, Simple2DCell::Ncells, access_mode::overwrite, and GPUArray< T >::resize().

◆ setVertexPositions()

void Simple2DCell::setVertexPositions ( vector< Dscalar2 >  newVertexPositions)
inherited

Set vertex positions according to a user-specified vector.

Does not update any other lists – it is the user's responsibility to maintain topology, etc, when using this function.

References ArrayHandle< T >::data, GPUArray< T >::getNumElements(), access_location::host, Simple2DCell::Nvertices, access_mode::overwrite, GPUArray< T >::resize(), and Simple2DCell::vertexPositions.

◆ setCellVelocitiesMaxwellBoltzmann()

Dscalar Simple2DCell::setCellVelocitiesMaxwellBoltzmann ( Dscalar  T)
inherited

Set velocities via a temperature. The return value is the total kinetic energy.

Set the cell velocities by drawing from a Maxwell-Boltzmann distribution, and then make sure there is no net momentum. The return value is the total kinetic energy

◆ setVertexVelocitiesMaxwellBoltzmann()

Dscalar Simple2DCell::setVertexVelocitiesMaxwellBoltzmann ( Dscalar  T)
inherited

Set velocities via a temperature for the vertex degrees of freedom.

Set the vertex velocities by drawing from a Maxwell-Boltzmann distribution, and then make sure there is no net momentum. The return value is the total kinetic energy.

◆ setModuliUniform()

void Simple2DCell::setModuliUniform ( Dscalar  newKA,
Dscalar  newKP 
)
inherited

set uniform moduli for all cells

set all cell K_A, K_P preferences to uniform values. PLEASE NOTE that as an optimization this data is not actually used at the moment, but the code could be trivially altered to use this

◆ setCellTypeUniform()

void Simple2DCell::setCellTypeUniform ( int  i)
inherited

◆ setCellType()

void Simple2DCell::setCellType ( vector< int > &  types)
inherited

Set cells to different "type".

Parameters
typesa vector of integers that the cell types will be set to

References Simple2DCell::cellType, ArrayHandle< T >::data, access_location::host, Simple2DCell::Ncells, access_mode::overwrite, and GPUArray< T >::resize().

◆ setVertexTopologyFromCells()

void Simple2DCell::setVertexTopologyFromCells ( vector< vector< int > >  cellVertexIndices)
inherited

An uncomfortable function to allow the user to set vertex topology "by hand".

This function allows a user to set the vertex topology by hand. The user is responsible for making sure the input topology is sensible. DMS NOTE – this functionality has not been thoroughly tested

Precondition
Nvertices and vertex positions are already set
Parameters
cellVertexIndicesa vector of vector of ints. Each vector of ints must correspond to the counter-clockwise ordering of vertices that make up the cell, and every vertex should appear at most three times in different cells

References Simple2DCell::cellVertexNum, Simple2DCell::cellVertices, ArrayHandle< T >::data, access_location::host, Simple2DCell::n_idx, Simple2DCell::Ncells, Simple2DCell::Nvertices, access_mode::overwrite, GPUArray< T >::resize(), Simple2DCell::vertexCellNeighbors, Simple2DCell::vertexMax, and Simple2DCell::vertexNeighbors.

◆ initializeCellSorting()

void Simple2DCell::initializeCellSorting ( )
protectedinherited

set the size of the cell-sorting structures, initialize lists simply

Sets the size of itt, tti, idxToTag, and tagToIdx, and sets all of them so that array[i] = i, i.e., unsorted

Precondition
Ncells is determined

References Simple2DCell::idxToTag, Simple2DCell::itt, Simple2DCell::Ncells, Simple2DCell::tagToIdx, and Simple2DCell::tti.

◆ initializeVertexSorting()

void Simple2DCell::initializeVertexSorting ( )
protectedinherited

set the size of the vertex-sorting structures, initialize lists simply

Sets the size of ittVertex, ttiVertex, idxToTagVertex, and tagToIdxVertex,and sets all of them so that array[i] = i, i.e., things are unsorted

Precondition
Nvertices is determined

References Simple2DCell::idxToTagVertex, Simple2DCell::ittVertex, Simple2DCell::Nvertices, Simple2DCell::tagToIdxVertex, and Simple2DCell::ttiVertex.

◆ reIndexCellArray() [1/3]

void Simple2DCell::reIndexCellArray ( GPUArray< int > &  array)
protectedinherited

◆ reIndexCellArray() [2/3]

void Simple2DCell::reIndexCellArray ( GPUArray< Dscalar > &  array)
protectedinherited

why use templates when you can type more?

Re-indexes GPUarrays of Dscalars

References ArrayHandle< T >::data, access_location::host, Simple2DCell::itt, Simple2DCell::Ncells, access_mode::read, and access_mode::readwrite.

◆ reIndexCellArray() [3/3]

void Simple2DCell::reIndexCellArray ( GPUArray< Dscalar2 > &  array)
protectedinherited

why use templates when you can type more?

Always called after spatial sorting is performed, reIndexCellArray shuffles the order of an array based on the spatial sort order of the cells

References ArrayHandle< T >::data, access_location::host, Simple2DCell::itt, Simple2DCell::Ncells, access_mode::read, and access_mode::readwrite.

◆ reIndexVertexArray()

void Simple2DCell::reIndexVertexArray ( GPUArray< Dscalar2 > &  array)
protectedinherited

why use templates when you can type more?

Called if the vertices need to be spatially sorted, reIndexVertexArray shuffles the order of an array based on the spatial sort order of the vertices

References ArrayHandle< T >::data, access_location::host, Simple2DCell::ittVertex, Simple2DCell::Nvertices, access_mode::read, and access_mode::readwrite.

◆ spatiallySortCells()

void Simple2DCell::spatiallySortCells ( )
protectedinherited

◆ spatiallySortVertices()

void Simple2DCell::spatiallySortVertices ( )
protectedinherited

◆ reportMeanCellForce()

void Simple2DCell::reportMeanCellForce ( bool  verbose)
virtualinherited

Report the current average force on each cell.

a utility/testing function...output the currently computed mean net force to screen.

Parameters
verboseif true also print out the force on each cell

Implements Simple2DModel.

References Simple2DCell::cellForces, Simple2DCell::cellPositions, access_location::host, and access_mode::read.

◆ reportq()

Dscalar Simple2DCell::reportq ( )
virtualinherited

Report the average value of p/sqrt(A) for the cells in the system.

Returns the mean value of the shape parameter:

Implements Simple2DModel.

◆ reportVarq()

Dscalar Simple2DCell::reportVarq ( )
inherited

Report the variance of p/sqrt(A) for the cells in the system.

Returns the variance of the shape parameter:

◆ reportVarAP()

Dscalar2 Simple2DCell::reportVarAP ( )
inherited

Report the variance of A and P for the cells in the system.

Returns the variance of the A and P for the system:

◆ reportMeanP()

Dscalar Simple2DCell::reportMeanP ( )
inherited

Report the mean value of the perimeter.

Returns the mean value of the perimeter

Member Data Documentation

◆ NeighIdxs

GPUArray<int2> voronoiModelBase::NeighIdxs
protected

An array that holds (particle, neighbor_number) info to avoid intra-warp divergence in GPU -based force calculations that might be used by child classes

Referenced by VoronoiQuadraticEnergy::computeVoronoiForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiSimpleTensionForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiTensionForceSetsGPU(), enforceTopology(), resizeAndReset(), and updateNeighIdxs().

◆ globalOnly

bool voronoiModelBase::globalOnly
protected

When true, the CPU branch will execute global retriangulations through CGAL on every time step.

When running on the CPU, should only global retriangulations be performed, or should local test-and-updates still be performed? Depending on parameters simulated, performance here can be quite difference, since the circumcircle test itself is CPU expensive

Referenced by setCPU(), and testTriangulationCPU().

◆ delSets

GPUArray<int2> voronoiModelBase::delSets
protected

delSet.data[n_idx(nn,i)] are the previous and next consecutive delaunay neighbors

These are orientationally ordered, of point i (for use in computing forces on GPU)

Referenced by VoronoiQuadraticEnergy::computeVoronoiForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiSimpleTensionForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiTensionForceSetsGPU(), enforceTopology(), getDelSets(), and resetLists().

◆ delOther

GPUArray<int> voronoiModelBase::delOther
protected

delOther.data[n_idx(nn,i)] contains the index of the "other" delaunay neighbor.

i.e., the mutual neighbor of delSet.data[n_idx(nn,i)].y and delSet.data[n_idx(nn,i)].z that isn't point i

Referenced by VoronoiQuadraticEnergy::computeVoronoiForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiSimpleTensionForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiTensionForceSetsGPU(), enforceTopology(), getDelSets(), and resetLists().

◆ vertexNeighbors

GPUArray<int> Simple2DCell::vertexNeighbors
inherited

VERTEX neighbors of every vertex.

in general, we have: vertexNeighbors[3*i], vertexNeighbors[3*i+1], and vertexNeighbors[3*i+2] contain the indices of the three vertices that are connected to vertex i

Referenced by vertexModelBase::cellDeath(), vertexModelBase::flipEdgesGPU(), Simple2DCell::setVertexTopologyFromCells(), Simple2DCell::spatiallySortVertices(), vertexModelBase::testAndPerformT1TransitionsCPU(), and vertexModelBase::testEdgesForT1GPU().

◆ vertexCellNeighbors

GPUArray<int> Simple2DCell::vertexCellNeighbors
inherited

◆ cellVertexNum

GPUArray<int> Simple2DCell::cellVertexNum
inherited

◆ cellType

GPUArray<int> Simple2DCell::cellType
inherited

An array of integers labeling cell type...an easy way of determining if cells are different.

Please note that "type" is not meaningful unless it is used by child classes. That is, things like area/perimeter preferences, or motility, or whatever are neither set nor accessed by cell type, but rather by cell index! Thus, this is just an additional data structure that can be useful. For instance, the VoronoiTension2D classes uses the integers of cellType to determine when to apply an additional line tension between cells.

Referenced by Simple2DCell::cellDeath(), Simple2DCell::cellDivision(), VertexQuadraticEnergyWithTension::computeVertexTensionForcesCPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiSimpleTensionForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiTensionForceSetsGPU(), Simple2DCell::setCellType(), Simple2DCell::setCellTypeUniform(), Simple2DCell::spatiallySortCells(), and Simple2DCell::spatiallySortVertices().

◆ tagToIdx

vector<int> Simple2DCell::tagToIdx
inherited

To write consistent files...the cell that started the simulation as index i has current index tagToIdx[i].

The Hilbert sorting stuff makes keeping track of particles, and re-indexing things when particle number changes, a pain. Here's a description of the four relevant data structures. tagToIdx[i] = a. At the beginning of a simulation, a particle had index "i", meaning its current state was found in position "i" of all various data vectors and arrays. That same particle's data is now in position "a" of those data structures. Short version: "Where do I look to find info for what I orinally called partice i?" idxToTag[a] = i. That is, idxToTag just helps invert the tagToIdx list. idxToTag[tagToIdx[i]]=i The above two structures (and the vertex versions of them) tell you how to go back and forth between the current state of the system and the initial state of the system. What about going back and forth between the current sorted state and the previous sorted state? The "itt" and "tti" vectors give this information. The itt and tti vectors are completely overwritten each time a spatial sorting is called. By the way, I apologize if the nomenclature of "index" vs. "tag" is the opposite of what you, the reader of these code comments, might expect.

Referenced by Simple2DCell::cellDeath(), Simple2DCell::cellDivision(), Simple2DCell::initializeCellSorting(), Simple2DCell::spatiallySortCells(), and Simple2DCell::spatiallySortVertices().

◆ cellVertices

GPUArray<int> Simple2DCell::cellVertices
protectedinherited

◆ voroCur

GPUArray<Dscalar2> Simple2DCell::voroCur
protectedinherited

3*Nvertices length array of the position of vertices around cells

For both vertex and Voronoi models, it may help to save the relative position of the vertices around a cell, either for easy force computation or in the geometry routine, etc. voroCur.data[n_idx(nn,i)] gives the nth vertex, in CCW order, of cell i

Referenced by VertexQuadraticEnergy::computeForcesCPU(), VertexQuadraticEnergy::computeForcesGPU(), vertexModelBase::computeGeometryCPU(), computeGeometryCPU(), vertexModelBase::computeGeometryGPU(), computeGeometryGPU(), VertexQuadraticEnergyWithTension::computeVertexTensionForcesCPU(), VoronoiQuadraticEnergy::computeVoronoiForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiSimpleTensionForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiTensionForceSetsGPU(), and resetLists().

◆ voroLastNext

GPUArray<Dscalar4> Simple2DCell::voroLastNext
protectedinherited

◆ itt

vector<int> Simple2DCell::itt
protectedinherited

A map between cell index and the spatially sorted version.

sortedArray[i] = unsortedArray[itt[i]] after a hilbert sort

Referenced by Simple2DCell::cellDeath(), Simple2DCell::cellDivision(), Simple2DCell::initializeCellSorting(), Simple2DCell::reIndexCellArray(), Simple2DCell::returnItt(), Simple2DCell::spatiallySortCells(), and Simple2DCell::spatiallySortVertices().


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