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

Implement a 2D active vertex model, using kernels in avmKernels. More...

#include <vertexQuadraticEnergy.h>

Inheritance diagram for VertexQuadraticEnergy:
Inheritance graph
[legend]

Public Member Functions

 VertexQuadraticEnergy (int n, Dscalar A0, Dscalar P0, bool reprod=false, bool runSPVToInitialize=false)
 the constructor: initialize as a Delaunay configuration with random positions and set all cells to have uniform target A_0 and P_0 parameters More...
 
virtual void computeForces ()
 compute the geometry and get the forces More...
 
virtual Dscalar computeEnergy ()
 compute the quadratic energy functional More...
 
void computeForcesCPU ()
 Compute the geometry (area & perimeter) of the cells on the CPU. More...
 
void computeForcesGPU ()
 Compute the geometry (area & perimeter) of the cells on the GPU. More...
 
virtual int getNumberOfDegreesOfFreedom ()
 In vertex models the number of degrees of freedom is the number of vertices.
 
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
 
void initializeVertexModelBase (int n, bool spvInitialize=false)
 Initialize vertexModelBase, set random orientations for vertex directors, prepare data structures. More...
 
void setCellsVoronoiTesselation (bool spvInitialize=false)
 Initialize cells to be a voronoi tesselation of a random point set. More...
 
virtual GPUArray< Dscalar2 > & returnForces ()
 return a reference to the GPUArray of the current forces
 
virtual GPUArray< Dscalar2 > & returnVelocities ()
 return a reference to the GPUArray of the current velocities
 
virtual GPUArray< Dscalar2 > & returnPositions ()
 return a reference to the GPUArray of the current positions
 
virtual GPUArray< Dscalar > & returnMasses ()
 return a reference to the GPUArray of the current masses
 
virtual void computeGeometryCPU ()
 Compute the geometry (area & perimeter) of the cells on the CPU. More...
 
virtual void computeGeometryGPU ()
 Compute the geometry (area & perimeter) of the cells on the GPU. More...
 
void getCellCentroids ()
 Call the CPU or GPU getCellCentroids function. More...
 
void getCellCentroidsCPU ()
 Get the cell position from the vertices on the CPU. More...
 
void getCellCentroidsGPU ()
 Get the cell position from the vertices on the GPU. More...
 
void getCellPositions ()
 Call the CPU or GPU getCellPositions function. More...
 
void getCellPositionsCPU ()
 Get the cell position from the average vertex position on the CPU. More...
 
void getCellPositionsGPU ()
 Get the cell position from the average vertex position on the GPU. More...
 
virtual void cellDivision (const vector< int > &parameters, const vector< Dscalar > &dParams={})
 Divide cell...vector should be cell index i, vertex 1 and vertex 2. More...
 
virtual void cellDeath (int cellIndex)
 Kill the indexed cell...cell must have only three associated vertices. More...
 
virtual void setT1Threshold (Dscalar t1t)
 Set the length threshold for T1 transitions.
 
void testAndPerformT1TransitionsCPU ()
 Simple test for T1 transitions (edge length less than threshold) on the CPU. More...
 
void testAndPerformT1TransitionsGPU ()
 Simple test for T1 transitions (edge length less than threshold) on the GPU...calls the following functions. More...
 
virtual void spatialSorting ()
 spatially sort the vertices along a Hilbert curve for data locality More...
 
virtual void enforceTopology ()
 update/enforce the topology, performing simple T1 transitions More...
 
void setCPU (bool global=true)
 Enforce CPU-only operation.
 
virtual void setCPU ()
 Enforce CPU-only operation. Derived classes might have to do more work when the CPU mode is invoked.
 
void reportNeighborsCell (int i)
 Handy for debugging T1 transitions...report the vertices owned by cell i.
 
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 computeGeometry ()
 call either the computeGeometryCPU or GPU routines for the current model More...
 
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< 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

GPUArray< int > vertexEdgeFlips
 flags that indicate whether an edge should be GPU-flipped (1) or not (0) More...
 
GPUArray< int > vertexEdgeFlipsCurrent
 it is important to not flip edges concurrently, so this data structure helps flip edges sequentially
 
GPUArray< Dscalar2 > vertexForceSets
 an array containing the three contributions to the force on each vertex More...
 
Dscalar T1Threshold
 A threshold defining the edge length below which a T1 transition will occur.
 
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 growCellVerticesList (int newVertexMax)
 if the maximum number of vertices per cell increases, grow the cellVertices list More...
 
void initializeEdgeFlipLists ()
 Initialize the data structures for edge flipping...should also be called if Nvertices changes. More...
 
void testEdgesForT1GPU ()
 test the edges for a T1 event, and grow the cell-vertex list if necessary More...
 
void flipEdgesGPU ()
 perform the edge flips found in the previous step More...
 
void getCellVertexSetForT1 (int v1, int v2, int4 &cellSet, int4 &vertexSet, bool &growList)
 For finding T1s on the CPU; find the set of vertices and cells involved in the transition. 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

GPUArray< int > growCellVertexListAssist
 data structure to help with cell-vertex list
 
GPUArray< int > finishedFlippingEdges
 data structure to help with not simultaneously trying to flip nearby edges
 
GPUArray< int > cellEdgeFlips
 data structure per cell for not simulataneously flipping nearby edges
 
GPUArray< int4 > cellSets
 data structure per cell for not simulataneously flipping nearby edges
 
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.
 

Friends

class AVMDatabaseNetCDF
 

Detailed Description

Implement a 2D active vertex model, using kernels in avmKernels.

A class that implements a simple active vertex model in 2D. This involves calculating forces on vertices, moving them around, and updating the topology of the cells according to some criteria.

This class is a child of the vertexModelBase class, which provides data structures like the positions of cells, vertex positions, indices of vertices around each cell, cells around each vertex, etc. updates/enforces the topology according to vertexModelBase' T1 functions

Constructor & Destructor Documentation

◆ VertexQuadraticEnergy()

VertexQuadraticEnergy::VertexQuadraticEnergy ( int  n,
Dscalar  A0,
Dscalar  P0,
bool  reprod = false,
bool  runSPVToInitialize = false 
)

the constructor: initialize as a Delaunay configuration with random positions and set all cells to have uniform target A_0 and P_0 parameters

Parameters
nnumber of CELLS to initialize
A0set uniform preferred area for all cells
P0set uniform preferred perimeter for all cells
reprodshould the simulation be reproducible (i.e. call a RNG with a fixed seed)
runSPVToInitializethe default constructor has the cells start as a Voronoi tesselation of a random point set. Set this flag to true to relax this initial configuration via the Voronoi2D class
Postcondition
Initialize(n,runSPVToInitialize) is called, setCellPreferencesUniform(A0,P0), and setModuliUniform(1.0,1.0)

Member Function Documentation

◆ computeForces()

void VertexQuadraticEnergy::computeForces ( )
virtual

compute the geometry and get the forces

compute the geometry and the forces and the vertices, on either the GPU or CPU as determined by flags

Reimplemented from Simple2DCell.

Reimplemented in VertexQuadraticEnergyWithTension.

References computeForcesCPU(), computeForcesGPU(), Simple2DCell::computeGeometry(), Simple2DCell::forcesUpToDate, and Simple2DCell::GPUcompute.

◆ computeEnergy()

Dscalar VertexQuadraticEnergy::computeEnergy ( )
virtual

compute the quadratic energy functional

Returns the quadratic energy functional: E = {cells} K_A(A_i-A_i,0)^2 + K_P(P_i-P_i,0)^2

Reimplemented from Simple2DCell.

Reimplemented in VertexQuadraticEnergyWithTension.

◆ computeForcesCPU()

void VertexQuadraticEnergy::computeForcesCPU ( )

Compute the geometry (area & perimeter) of the cells on the CPU.

Use the data pre-computed in the geometry routine to rapidly compute the net force on each vertex

References Simple2DCell::AreaPeri, Simple2DCell::AreaPeriPreferences, access_location::host, access_mode::overwrite, access_mode::read, Simple2DCell::vertexCellNeighbors, Simple2DCell::vertexForces, vertexModelBase::vertexForceSets, Simple2DCell::voroCur, and Simple2DCell::voroLastNext.

Referenced by VertexQuadraticEnergyWithTension::computeForces(), and computeForces().

◆ computeForcesGPU()

void VertexQuadraticEnergy::computeForcesGPU ( )

◆ moveDegreesOfFreedom()

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

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

move vertices according to an inpute GPUarray

Reimplemented from Simple2DCell.

Referenced by vertexModelBase::getNumberOfDegreesOfFreedom().

◆ initializeVertexModelBase()

void vertexModelBase::initializeVertexModelBase ( int  n,
bool  spvInitialize = false 
)
inherited

Initialize vertexModelBase, set random orientations for vertex directors, prepare data structures.

Take care of all base class initialization functions, this involves setting arrays to the right size, etc.

Referenced by vertexModelBase::getForces().

◆ setCellsVoronoiTesselation()

void vertexModelBase::setCellsVoronoiTesselation ( bool  spvInitialize = false)
inherited

Initialize cells to be a voronoi tesselation of a random point set.

A function of convenience.... initialize cell positions and vertices by constructing the Delaunay triangulation of the current cell positions. If you want something more regular, run the Voronoi mode for a few timesteps to smooth out the random point set first.

Parameters
spvInitializeonly use if the initial cell positions are to be random, and you want to make the points more uniform
Postcondition
After this is called, all topology data structures are initialized

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

◆ computeGeometryCPU()

void vertexModelBase::computeGeometryCPU ( )
virtualinherited

Compute the geometry (area & perimeter) of the cells on the CPU.

Very similar to the function in Voronoi2d.cpp, but optimized since we already have some data structures (the vertices)...compute the area and perimeter of the cells

Reimplemented from Simple2DCell.

References Simple2DCell::AreaPeri, Simple2DCell::cellVertexNum, Simple2DCell::cellVertices, ArrayHandle< T >::data, access_location::host, Simple2DCell::Ncells, access_mode::read, access_mode::readwrite, Simple2DCell::vertexCellNeighbors, Simple2DCell::vertexPositions, Simple2DCell::voroCur, and Simple2DCell::voroLastNext.

Referenced by vertexModelBase::returnMasses().

◆ computeGeometryGPU()

void vertexModelBase::computeGeometryGPU ( )
virtualinherited

Compute the geometry (area & perimeter) of the cells on the GPU.

Very similar to the function in Voronoi2d.cpp, but optimized since we already have some data structures (the vertices)

Reimplemented from Simple2DCell.

References Simple2DCell::AreaPeri, Simple2DCell::Box, Simple2DCell::cellVertexNum, Simple2DCell::cellVertices, ArrayHandle< T >::data, access_location::device, gpu_vm_geometry(), Simple2DCell::n_idx, Simple2DCell::Ncells, access_mode::overwrite, access_mode::read, Simple2DCell::vertexCellNeighbors, Simple2DCell::vertexPositions, Simple2DCell::voroCur, and Simple2DCell::voroLastNext.

Referenced by vertexModelBase::returnMasses().

◆ getCellCentroids()

void vertexModelBase::getCellCentroids ( )
inherited

Call the CPU or GPU getCellCentroids function.

This function fills the "cellPositions" GPUArray with the centroid of every cell. Does not assume that the area in the AreaPeri array is current. This function just calls the CPU or GPU routine, as determined by the GPUcompute flag

References vertexModelBase::getCellCentroidsCPU(), vertexModelBase::getCellCentroidsGPU(), and Simple2DCell::GPUcompute.

Referenced by vertexModelBase::returnMasses().

◆ getCellCentroidsCPU()

void vertexModelBase::getCellCentroidsCPU ( )
inherited

◆ getCellCentroidsGPU()

void vertexModelBase::getCellCentroidsGPU ( )
inherited

Get the cell position from the vertices on the GPU.

GPU computation of the centroid of every cell

Referenced by vertexModelBase::getCellCentroids(), and vertexModelBase::returnMasses().

◆ getCellPositions()

void vertexModelBase::getCellPositions ( )
inherited

Call the CPU or GPU getCellPositions function.

This function fills the "cellPositions" GPUArray with the mean position of the vertices of each cell. This function just calls the CPU or GPU routine, as determined by the GPUcompute flag

References vertexModelBase::getCellPositionsCPU(), vertexModelBase::getCellPositionsGPU(), and Simple2DCell::GPUcompute.

Referenced by vertexModelBase::returnMasses().

◆ getCellPositionsCPU()

void vertexModelBase::getCellPositionsCPU ( )
inherited

Get the cell position from the average vertex position on the CPU.

One would prefer the cell position to be defined as the centroid, requiring an additional computation of the cell area. This may be implemented some day, but for now we define the cell position as the straight average of the vertex positions. This isn't really used much, anyway, so update this only when the functionality becomes needed

References Simple2DCell::cellPositions, Simple2DCell::cellVertexNum, Simple2DCell::cellVertices, access_location::host, access_mode::read, access_mode::readwrite, and Simple2DCell::vertexPositions.

Referenced by vertexModelBase::getCellPositions(), and vertexModelBase::returnMasses().

◆ getCellPositionsGPU()

void vertexModelBase::getCellPositionsGPU ( )
inherited

◆ cellDivision()

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

Divide cell...vector should be cell index i, vertex 1 and vertex 2.

Trigger a cell division event, which involves some laborious re-indexing of various data structures. This simple version of cell division will take a cell and two specified vertices. The edges emanating clockwise from each of the two vertices will gain a new vertex in the middle of those edges. A new cell is formed by connecting those two new vertices together. The vector of "parameters" here should be three integers: parameters[0] = the index of the cell to undergo a division event parameters[1] = the first vertex to gain a new (clockwise) vertex neighbor. parameters[2] = the second ..... The two vertex numbers should be between 0 and celLVertexNum[parameters[0]], respectively, NOT the indices of the vertices being targeted Note that dParams does nothing

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

Reimplemented from Simple2DActiveCell.

References Simple2DCell::Ncells.

Referenced by vertexModelBase::returnMasses().

◆ cellDeath()

void vertexModelBase::cellDeath ( int  cellIndex)
virtualinherited

Kill the indexed cell...cell must have only three associated vertices.

Trigger a cell death event. This REQUIRES that the vertex model cell to die be a triangle (i.e., we are mimicking a T2 transition)

Reimplemented from Simple2DActiveCell.

References Simple2DCell::cellVertexNum, Simple2DCell::cellVertices, ArrayHandle< T >::data, Simple2DCell::n_idx, Simple2DCell::Nvertices, Simple2DCell::vertexCellNeighbors, and Simple2DCell::vertexNeighbors.

Referenced by vertexModelBase::returnMasses().

◆ testAndPerformT1TransitionsCPU()

void vertexModelBase::testAndPerformT1TransitionsCPU ( )
inherited

Simple test for T1 transitions (edge length less than threshold) on the CPU.

Test whether a T1 needs to be performed on any edge by simply checking if the edge length is beneath a threshold. This function also performs the transition and maintains the auxiliary data structures

References Simple2DCell::cellVertexNum, Simple2DCell::cellVertices, access_location::host, access_mode::readwrite, Simple2DCell::vertexCellNeighbors, Simple2DCell::vertexNeighbors, and Simple2DCell::vertexPositions.

Referenced by vertexModelBase::enforceTopology().

◆ testAndPerformT1TransitionsGPU()

void vertexModelBase::testAndPerformT1TransitionsGPU ( )
inherited

Simple test for T1 transitions (edge length less than threshold) on the GPU...calls the following functions.

Because the cellVertexList might need to grow, it's convenient to break this into two parts

References vertexModelBase::flipEdgesGPU(), and vertexModelBase::testEdgesForT1GPU().

Referenced by vertexModelBase::enforceTopology().

◆ spatialSorting()

void vertexModelBase::spatialSorting ( )
virtualinherited

spatially sort the vertices along a Hilbert curve for data locality

When sortPeriod < 0 this routine does not get called

Postcondition
vertices are re-ordered according to a Hilbert sorting scheme, cells are reordered according to what vertices they are near, and all data structures are updated

Reimplemented from Simple2DModel.

References Simple2DCell::reIndexVertexArray(), Simple2DActiveCell::spatiallySortVerticesAndCellActivity(), Simple2DCell::vertexMasses, and Simple2DCell::vertexVelocities.

◆ enforceTopology()

void vertexModelBase::enforceTopology ( )
virtualinherited

update/enforce the topology, performing simple T1 transitions

enforce and update topology of vertex wiring on either the GPU or CPU

Reimplemented from Simple2DCell.

References Simple2DCell::GPUcompute, vertexModelBase::testAndPerformT1TransitionsCPU(), and vertexModelBase::testAndPerformT1TransitionsGPU().

◆ growCellVerticesList()

void vertexModelBase::growCellVerticesList ( int  newVertexMax)
protectedinherited

if the maximum number of vertices per cell increases, grow the cellVertices list

when a transition increases the maximum number of vertices around any cell in the system, call this function first to copy over the cellVertices structure into a larger array

References Simple2DCell::cellVertexNum, Simple2DCell::cellVertices, ArrayHandle< T >::data, access_location::host, Simple2DCell::n_idx, Simple2DCell::Ncells, access_mode::read, access_mode::readwrite, GPUArray< T >::resize(), and Simple2DCell::vertexMax.

Referenced by vertexModelBase::setCPU(), and vertexModelBase::testEdgesForT1GPU().

◆ initializeEdgeFlipLists()

void vertexModelBase::initializeEdgeFlipLists ( )
protectedinherited

Initialize the data structures for edge flipping...should also be called if Nvertices changes.

Initialize the auxilliary edge flip data structures to zero

References ArrayHandle< T >::data, vertexModelBase::finishedFlippingEdges, access_location::host, Simple2DCell::Nvertices, access_mode::overwrite, GPUArray< T >::resize(), vertexModelBase::vertexEdgeFlips, and vertexModelBase::vertexEdgeFlipsCurrent.

Referenced by vertexModelBase::setCPU().

◆ testEdgesForT1GPU()

void vertexModelBase::testEdgesForT1GPU ( )
protectedinherited

◆ flipEdgesGPU()

void vertexModelBase::flipEdgesGPU ( )
protectedinherited

◆ getCellVertexSetForT1()

void vertexModelBase::getCellVertexSetForT1 ( int  vertex1,
int  vertex2,
int4 &  cellSet,
int4 &  vertexSet,
bool &  growList 
)
protectedinherited

For finding T1s on the CPU; find the set of vertices and cells involved in the transition.

A utility function for the CPU T1 transition routine. Given two vertex indices representing an edge that will undergo a T1 transition, return in the pass-by-reference variables a helpful representation of the cells in the T1 and the vertices to be re-wired...see the comments in "testAndPerformT1TransitionsCPU" for what that representation is

References Simple2DCell::cellVertexNum, Simple2DCell::cellVertices, ArrayHandle< T >::data, access_location::host, Simple2DCell::n_idx, access_mode::read, Simple2DCell::vertexCellNeighbors, and Simple2DCell::vertexMax.

◆ 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 voronoiModelBase::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

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 Simple2DCell::computeGeometryCPU(), Simple2DCell::computeGeometryGPU(), and Simple2DCell::GPUcompute.

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

◆ 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

◆ vertexEdgeFlips

GPUArray<int> vertexModelBase::vertexEdgeFlips
inherited

flags that indicate whether an edge should be GPU-flipped (1) or not (0)

if vertexEdgeFlips[3*i+j]=1 (where j runs from 0 to 2), the the edge connecting vertex i and vertex vertexNeighbors[3*i+j] has been marked for a T1 transition

Referenced by vertexModelBase::flipEdgesGPU(), vertexModelBase::initializeEdgeFlipLists(), and vertexModelBase::testEdgesForT1GPU().

◆ vertexForceSets

GPUArray<Dscalar2> vertexModelBase::vertexForceSets
inherited

an array containing the three contributions to the force on each vertex

vertexForceSets[3*i], vertexForceSets[3*i+1], and vertexForceSets[3*i+2] contain the contribution to the net force on vertex i due to the three cell neighbors of vertex i

Referenced by computeForcesCPU(), computeForcesGPU(), and VertexQuadraticEnergyWithTension::computeVertexTensionForcesCPU().

◆ 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 computeForcesCPU(), 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
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: