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

Implement data structures and functions common to many off-lattice models of cells in 2D. More...

#include <Simple2DCell.h>

Inheritance diagram for Simple2DCell:
Inheritance graph
[legend]

Public Member Functions

 Simple2DCell ()
 initialize member variables to some defaults More...
 
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 int getNumberOfDegreesOfFreedom ()
 get the number of degrees of freedom, defaulting to the number of cells
 
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 void computeGeometryCPU ()
 let computeGeometryCPU be defined in derived classes
 
virtual void computeGeometryGPU ()
 let computeGeometryGPU be defined in derived classes
 
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...
 
virtual void getForces (GPUArray< Dscalar2 > &forces)
 copy the models current set of forces to the variable
 
virtual void moveDegreesOfFreedom (GPUArray< Dscalar2 > &displacements, Dscalar scale=1.)
 move the degrees of freedom
 
virtual void enforceTopology ()
 Do everything necessary to update or enforce the topology in the current model.
 
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...
 
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)
 allow for cell death, killing off the cell with the specified index 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< Dscalar2 > & returnForces ()
 Return a reference to forces 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 setCPU (bool a)=0
 
virtual void setv0Dr (Dscalar a, Dscalar b)=0
 
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 spatialSorting ()
 do everything necessary to perform a Hilbert sort
 
virtual void setTime (Dscalar time)
 set the time
 

Public Attributes

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 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

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

Implement data structures and functions common to many off-lattice models of cells in 2D.

A class defining some of the fundamental attributes and operations common to 2D off-lattice models of cells. Note that while all 2D off-lattice models use some aspects of this base class, not all of them are required to implement or use all of the below

Constructor & Destructor Documentation

◆ Simple2DCell()

Simple2DCell::Simple2DCell ( )

initialize member variables to some defaults

An extremely simple constructor that does nothing, but enforces default GPU operation

References Box, and forcesUpToDate.

Member Function Documentation

◆ initializeSimple2DCell()

void Simple2DCell::initializeSimple2DCell ( int  n)

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 ( )
virtual

call either the computeGeometryCPU or GPU routines for the current model

Simply call either the CPU or GPU routine in the current or derived model

References computeGeometryCPU(), computeGeometryGPU(), and GPUcompute.

Referenced by VertexQuadraticEnergyWithTension::computeForces(), VoronoiQuadraticEnergyWithTension::computeForces(), VertexQuadraticEnergy::computeForces(), VoronoiQuadraticEnergy::computeForces(), and computeForces().

◆ computeKineticEnergy()

Dscalar Simple2DCell::computeKineticEnergy ( )

Call masses and velocities to get the total kinetic energy.

E = 0.5*m_i v_i^2

◆ computeKineticPressure()

Dscalar4 Simple2DCell::computeKineticPressure ( )

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 
)

Set uniform cell area and perimeter preferences.

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

Referenced by enforceTopology().

◆ setCellPreferences()

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

Set cell area and perimeter preferences according to input vector.

Set the Area and Perimeter preferences to the input vector

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

◆ setCellPositionsRandomly()

void Simple2DCell::setCellPositionsRandomly ( )

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 cellPositions, Ncells, and GPUArray< T >::resize().

◆ cellDivision()

void Simple2DCell::cellDivision ( const vector< int > &  parameters,
const vector< Dscalar > &  dParams = {} 
)
virtual

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

This function supports cellDivisions, updating data structures in Simple2DCell This function will grow the cell lists by 1 and assign the new cell (the last element of those arrays) the values of the cell given by parameters[0] Note that dParams does nothing by default, but allows more general virtual functions to be defined downstream (used in the Voronoi branch)

Reimplemented in vertexModelBase, voronoiModelBase, and Simple2DActiveCell.

References AreaPeri, AreaPeriPreferences, cellMasses, cellPositions, cellType, cellVelocities, ArrayHandle< T >::data, forcesUpToDate, idxToTag, itt, make_Dscalar2(), Moduli, n_idx, Ncells, GPUArray< T >::resize(), tagToIdx, tti, and vertexMax.

Referenced by Simple2DActiveCell::cellDivision().

◆ cellDeath()

void Simple2DCell::cellDeath ( int  cellIndex)
virtual

allow for cell death, killing off the cell with the specified index

When beggars die, there are no comets seen; The heavens themselves blaze forth the death of princes... Which are your cells? This function supports removing a single cell from the simulation, which requires re-indexing and relabeling the data structures in the simulation.

Postcondition
Simple2DCell data vectors are one element shorter, and cells with a higher index than "cellIndex" get re-labeled down by one.

Reimplemented in vertexModelBase, voronoiModelBase, and Simple2DActiveCell.

References AreaPeri, AreaPeriPreferences, cellMasses, cellPositions, cellType, cellVelocities, forcesUpToDate, idxToTag, itt, Moduli, Ncells, GPUArray< T >::resize(), tagToIdx, and tti.

Referenced by Simple2DActiveCell::cellDeath().

◆ setCellPositions()

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

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 cellPositions, ArrayHandle< T >::data, GPUArray< T >::getNumElements(), access_location::host, Ncells, access_mode::overwrite, and GPUArray< T >::resize().

◆ setVertexPositions()

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

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, Nvertices, access_mode::overwrite, GPUArray< T >::resize(), and vertexPositions.

◆ setCellVelocitiesMaxwellBoltzmann()

Dscalar Simple2DCell::setCellVelocitiesMaxwellBoltzmann ( Dscalar  T)

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)

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 
)

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)

Set all cells to the same "type".

set all cell types to i

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

◆ setCellType()

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

Set cells to different "type".

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

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

◆ setVertexTopologyFromCells()

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

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 cellVertexNum, cellVertices, ArrayHandle< T >::data, access_location::host, n_idx, Ncells, Nvertices, access_mode::overwrite, GPUArray< T >::resize(), vertexCellNeighbors, vertexMax, and vertexNeighbors.

◆ initializeCellSorting()

void Simple2DCell::initializeCellSorting ( )
protected

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 idxToTag, itt, Ncells, tagToIdx, and tti.

◆ initializeVertexSorting()

void Simple2DCell::initializeVertexSorting ( )
protected

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 idxToTagVertex, ittVertex, Nvertices, tagToIdxVertex, and ttiVertex.

◆ reIndexCellArray() [1/3]

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

◆ reIndexCellArray() [2/3]

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

why use templates when you can type more?

Re-indexes GPUarrays of Dscalars

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

◆ reIndexCellArray() [3/3]

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

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, itt, Ncells, access_mode::read, and access_mode::readwrite.

◆ reIndexVertexArray()

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

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, ittVertex, Nvertices, access_mode::read, and access_mode::readwrite.

◆ spatiallySortCells()

void Simple2DCell::spatiallySortCells ( )
protected

Perform a spatial sorting of the cells to try to maintain data locality.

take the current location of the cells and sort them according the their order along a 2D Hilbert curve

References AreaPeri, AreaPeriPreferences, Box, cellMasses, cellPositions, cellType, cellVelocities, ArrayHandle< T >::data, HilbertSorter::getIdx(), access_location::host, idxToTag, itt, Moduli, Ncells, access_mode::readwrite, reIndexCellArray(), tagToIdx, and tti.

Referenced by Simple2DActiveCell::spatiallySortCellsAndCellActivity().

◆ spatiallySortVertices()

void Simple2DCell::spatiallySortVertices ( )
protected

Perform a spatial sorting of the vertices to try to maintain data locality.

take the current location of the vertices and sort them according the their order along a 2D Hilbert curve. This routine first sorts the vertices, and then uses the vertex sorting to derive a sorting of the cells

Postcondition
both the itt, tti,... and ittVertex, ttiVertex... arrays are correctly set

References AreaPeriPreferences, Box, cellPositions, cellType, cellVertexNum, cellVertices, ArrayHandle< T >::data, HilbertSorter::getIdx(), access_location::host, idxToTag, idxToTagVertex, itt, ittVertex, Moduli, n_idx, Ncells, Nvertices, access_mode::read, access_mode::readwrite, reIndexCellArray(), reIndexVertexArray(), tagToIdx, tagToIdxVertex, tti, ttiVertex, vertexCellNeighbors, vertexMasses, vertexNeighbors, vertexPositions, and vertexVelocities.

Referenced by Simple2DActiveCell::spatiallySortVerticesAndCellActivity().

◆ reportMeanCellForce()

void Simple2DCell::reportMeanCellForce ( bool  verbose)
virtual

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 cellForces, cellPositions, access_location::host, and access_mode::read.

◆ reportq()

Dscalar Simple2DCell::reportq ( )
virtual

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 ( )

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

Returns the variance of the shape parameter:

◆ reportVarAP()

Dscalar2 Simple2DCell::reportVarAP ( )

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 ( )

Report the mean value of the perimeter.

Returns the mean value of the perimeter

Member Data Documentation

◆ vertexNeighbors

GPUArray<int> Simple2DCell::vertexNeighbors

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(), setVertexTopologyFromCells(), spatiallySortVertices(), vertexModelBase::testAndPerformT1TransitionsCPU(), and vertexModelBase::testEdgesForT1GPU().

◆ vertexCellNeighbors

GPUArray<int> Simple2DCell::vertexCellNeighbors

◆ cellVertexNum

GPUArray<int> Simple2DCell::cellVertexNum

◆ cellType

GPUArray<int> Simple2DCell::cellType

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 cellDeath(), cellDivision(), VertexQuadraticEnergyWithTension::computeVertexTensionForcesCPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiSimpleTensionForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiTensionForceSetsGPU(), setCellType(), setCellTypeUniform(), spatiallySortCells(), and spatiallySortVertices().

◆ tagToIdx

vector<int> Simple2DCell::tagToIdx

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 cellDeath(), cellDivision(), initializeCellSorting(), spatiallySortCells(), and spatiallySortVertices().

◆ cellVertices

GPUArray<int> Simple2DCell::cellVertices
protected

A structure that indexes the vertices defining each cell.

cellVertices is a large, 1D array containing the vertices associated with each cell. It must be accessed with the help of the Index2D structure n_idx. the index of the kth vertex of cell c (where the ordering is counter-clockwise starting with a random vertex) is given by cellVertices[n_idx(k,c)];

Referenced by vertexModelBase::cellDeath(), vertexModelBase::computeGeometryCPU(), vertexModelBase::computeGeometryGPU(), VertexQuadraticEnergyWithTension::computeVertexTensionForcesCPU(), vertexModelBase::flipEdgesGPU(), vertexModelBase::getCellCentroidsCPU(), vertexModelBase::getCellPositionsCPU(), vertexModelBase::getCellPositionsGPU(), vertexModelBase::getCellVertexSetForT1(), vertexModelBase::growCellVerticesList(), vertexModelBase::reportNeighborsCell(), setVertexTopologyFromCells(), spatiallySortVertices(), vertexModelBase::testAndPerformT1TransitionsCPU(), and vertexModelBase::testEdgesForT1GPU().

◆ voroCur

GPUArray<Dscalar2> Simple2DCell::voroCur
protected

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(), voronoiModelBase::computeGeometryCPU(), vertexModelBase::computeGeometryGPU(), voronoiModelBase::computeGeometryGPU(), VertexQuadraticEnergyWithTension::computeVertexTensionForcesCPU(), VoronoiQuadraticEnergy::computeVoronoiForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiSimpleTensionForceSetsGPU(), VoronoiQuadraticEnergyWithTension::computeVoronoiTensionForceSetsGPU(), and voronoiModelBase::resetLists().

◆ voroLastNext

GPUArray<Dscalar4> Simple2DCell::voroLastNext
protected

◆ itt

vector<int> Simple2DCell::itt
protected

A map between cell index and the spatially sorted version.

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

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


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