Source code documentation¶
General CPM and PDE code¶
-
class Dish¶
The virtual Petri dish. Hosts the cells with states and the CA-plane.
Public Functions
-
void Init(void)¶
Init defines the initial state of the virtual cell culture.
Define Init() in your main file describing the simulation set up, within the block INIT { }. See for examples vessel.cpp and sorting.cpp.
-
void Erase(void)¶
Erase all cells.
-
int Time(void) const¶
Returns the number of completed Monte Carlo Steps.
-
int CountCells(void) const¶
Returns the number of cells in the dish, excluding apoptosed cells.
-
void CellGrowthAndDivision(void)¶
Stretched induced cell growth and division.
See Hogeweg (2000), Journal of Theoretical Biology.
Find stretched cells, and increase their target area. Find enlarged cells, and divide them.
-
int Area(void) const¶
. Returns the summed area of all cells in the dish
-
int TargetArea(void) const¶
Returns the summed of all cells target area in the dish.
-
int SizeX(void)¶
Returns the horizontal size of the dish.
-
int SizeY(void)¶
Returns the horizontal size of the dish.
-
void Init(void)¶
-
class CellularPotts¶
Public Functions
-
int GetMatrixLevel(int x, int y)¶
Obtain the level of matrix.
- Returns:
Matrix concentration concentration
-
int GetActLevel(int x, int y)¶
Obtain the level of act.
- Returns:
Act concentration
-
CellularPotts(std::vector<Cell> *cells, const int sizex = 200, const int sizey = 200)¶
Constructs a CA field. This should be done in “Dish”.
-
void InitialiseEdgeList(void)¶
Initialises the edgelist at the beginning of a simulation.
The edgelist keeps track of pairs of lattice points that are eligible to change the CPM configuration. This function initialises the edgelist at the start.
-
virtual void AllocateSigma(int sx, int sy)¶
Allocates data for the sigma array.
Keyword virtual means, that derived classed (cppvmCellularPotts) can override this function and carry out the memory allocation in their preferred way Every time AllocateSigma is called in the base class methods the function belonging the actual type will be called
-
int **SearchNandPlot(Graphics *g = 0, bool get_neighbours = true)¶
Plots the dish to the screen or to a movie and searches the neighbours.
These distinct tasks have been lumped together in the same method because both for drawing the black lines between the cells and for searching the neighbours the cell borders have to be determined.
- Returns:
neighbourhood array
-
inline void Plot(Graphics *g)¶
Plot the dish to Graphics window g.
-
void PlotIsing(Graphics *g, int mag)¶
Special plotting for Ising model Only plot lines between the two states.
Plot in black & white for the Ising model
-
void SearchNandPlotClear(Graphics *g = 0)¶
Plot the neighbours between cells.
-
int **SearchNeighboursMatrix()¶
Obtain the neighbour matrix of cells.
- Returns:
Neighbourhood array
-
int GetNewPerimeterIfXYWereAdded(int sxyp, int x, int y)¶
Get perimeter if new pixel is added.
- Returns:
New perimeter
-
int GetNewPerimeterIfXYWereRemoved(int sxy, int x, int y)¶
Get perimeter if pixel is removed.
- Returns:
New perimeter
-
void SetTargetPerimeter(int tau, int value)¶
Set the target perimeter.
-
void SetLambdaPerimeter(int tau, int value)¶
Set the strength of the perimeter constraint.
-
inline int **SearchNeighbours(void)¶
Searches the cells’ neighbours without plotting.
-
inline int Mass(void)¶
Return the total area occupied by the cells.
-
void FindBoundingBox(void)¶
Find a bounding box that contains all cells Currently has no output.
-
void PlotSigma(Graphics *g, int mag = 2)¶
Plot the cells according to their cell identity, not their type. The black lines are omitted.
A simple method to plot all sigma’s in window without the black lines
-
inline void DivideCells(void)¶
Divide all cells. Divide along cell elongation axis.
-
void DivideCells(std::vector<bool> which_cells)¶
Divide all cells marked “true” in which_cells.
If which_cells is empty, this method divides all cells.
- Parameters:
which_cells – is a vector<bool> with the same number of elements as the number of cells. It is a mask indicating which cells should be divided; each cell marked true will be divided.
-
int AmoebaeMove(PDE *PDEfield = 0, bool anneal = false)¶
Monte Carlo Step. Returns summed energy change.
Implements the core CPM algorithm. Carries out one MCS.
- Returns:
Total energy change during MCS.
-
int Act_AmoebaeMove(PDE *PDEfield)¶
Implements the core CPM algorithm including Act dynamics. Carries out one MCS with the edge lsit algorithmAMo.
- Returns:
Total energy change during MCS.
-
int KawasakiMove(PDE *PDEfield = 0)¶
Monte Carlo Step. Returns summed energy change.
Implements the core CPM algorithm with Kawasaki dynamics. Carries out one MCS.
- Returns:
Total energy change during MCS.
-
int IsingMove(PDE *PDEfield = 0)¶
Monte Carlo Step. Returns summed energy change.
Implements Metropolis dynamics for the Ising model. Carries out one MCS.
- Returns:
Total energy change during MCS.
-
int PottsMove(PDE *PDEfield = 0)¶
Monte Carlo Step. Returns summed energy change.
Implements standard large q-Potts model. Carries out one MCS.
- Returns:
Total energy change during MCS.
-
int PottsNeighbourMove(PDE *PDEfield)¶
Monte Carlo Step. Returns summed energy change.
Implements standard large q-Potts model via Neighbour copies. Carries out one MCS.
- Returns:
Total energy change during MCS.
-
CellECMInteractions GetCellECMInteractions() const¶
Returns changes made to the adhesions since the last reset.
- Returns:
The accumulated changes
-
void ResetCellECMInteractions()¶
Clears recorded changes to the adhesions.
-
void SetECMBoundaryState(ECMBoundaryState const &ecm_boundary_state)¶
Set ECM boundary state, overwriting the current state.
-
void ReadZygotePicture(void)¶
Read initial cell shape from XPM file. Reads the initial cell shape from an include xpm picture called “ZYGXPM(ZYGOTE)”, and it allocates enough cells for it to the Dish.
-
void ConstructInitCells(Dish &beast)¶
Construct initial cells Construct the cells from the sigma-field.
-
inline int Time() const¶
Returns the number of completed Monte Carlo steps.
-
inline int SizeX() const¶
Return the horizontal size of the CA plane.
-
inline int SizeY() const¶
Return the vertical size of the CA plane.
-
inline int Sigma(const int x, const int y) const¶
Return the value of lattice site (x,y).
i.e. This will return the index of the cell which occupies site (x,y).
-
Dir *FindCellDirections(void) const¶
In this method the principal axes of the cells are computed using the method described in “Biometry”, box 15.5
- Returns:
a pointer to a “new[]”ed array containing the directions. The memory has to be freed afterwards using the delete[] operator
-
int ThrowInCells(int n, int cellsize)¶
Initialise the CA plane with n circular cells fitting in a cellsize^2 square.
! Fill the plane with initial cells
- Returns:
actual amount of cells (some are not draw due to overlap)
-
int GrowInCells(int n_cells, int cellsize, double subfield = 1., int posx = -1, int posy = -1)¶
Initialise the subfield in which eden growth takes place.
- Parameters:
n – Number of cells.
cellsize – Number of Eden growth iterations.
subfield – Defines a centered frame of size (size/subfield)^2 in which all cell will be positioned.
- Returns:
Index of last cell inserted.
-
int GrowInCells(int n_cells, int cell_size, int sx, int sy, int offset_x, int offset_y)¶
Initialise the CA plane with n cells using an Eden growth algorithm.
- Parameters:
n – Number of cells.
cell_size – Number of Eden growth iterations.
sx – x-size of subfield.
sy – y-size of subfield.
offset_x – x location for subfield.
offset_y – y location for subfield.
- Returns:
Index of last cell inserted.
-
void RandomSpins(double prob)¶
Initialise cpm field with a random sigma of 0 or 1 of every pixel.
- Parameters:
prob – This fraction of pixels will be a medium pixel
-
int SquareCell(int sig, int cx, int cy, int size)¶
Initialise a square cell.
Draw a square cell in at (cx,cy)
- Parameters:
sig – sigma value of the cell
cx – x-coordinate of the cell
cy – y-coordinate of the cell
size – length of the square cells
-
void ShowDirections(Graphics &g, const Dir *celldir) const¶
Display the division planes returned by FindCellDirections.
- Parameters:
g – Graphics window
celldir – cell axes as returned by FindCellDirections.
-
double MeanCellArea(void) const¶
Returns the mean area of the cells.
-
double CellDensity(void) const¶
Returns the cell density.
Cell density is defined as the area occupied by cells divided by the size of the field.
-
void ResetTargetLengths(void)¶
Set target lengths of all cells to the value given in parameter file.
-
void SetRandomTypes(void)¶
Give each cell a random cell type. The number of cell types is defined by the J parameter file. (See Jtable in parameter file).
-
void GrowAndDivideCells(int growth_rate)¶
Cells grow until twice their original target_length, then divide, with rate “growth_rate”
-
double DrawConvexHull(Graphics *g, int color = 1)¶
Draw convex hull around all cells.
- Returns:
The area of the convex hull in lattice sites.
-
double Compactness(void)¶
Calculate compactness (summed_area/hull_area) of all cells. This is a good measure for the density. Function allows a bounding box.
- Returns:
Compactness.
-
void RandomSigma(int n_cells)¶
Assign random sigma to every lattice point.
- n_cells: total number of cells
-
void MeasureCellSizes(void)¶
Measure the initial cell sizes Measure cell sizes of all initial size and assign them to the cells.
-
void MeasureCellPerimeters()¶
Measure the initial cell perimeters Measure cell perimeters of all initial size and assign them to the cells.
-
void anneal(int steps)¶
Run amoebaemove while only accepting negative delta H.
-
int **get_annealed_sigma(int steps)¶
Find the sigma field after annealing steps.
- Parameters:
steps – Number of annealing MCS
- Returns:
sigma-field after annealing
-
bool plotPos(int x, int y, Graphics *graphics)¶
plot the sigma at (x,y)
- Returns:
True if cell belongs to medium
-
void linePlotPos(int x, int y, Graphics *graphics)¶
plot cell outlines
-
int GetMatrixLevel(int x, int y)¶
-
class Cell¶
Public Functions
-
inline Cell(const Dish &who, int settau = 1)¶
Constructor to insert a cell into Dish “who” Used to add a new Cell to the dish: new Cell(dish,celtype).
-
void CellBirth(Cell &mother)¶
Add a new cell to the dish.
Call it as: new Cell(parent, true); mother will be modified for ancestry administration!
- Parameters:
settau – Cell type of daughter cell.
-
inline Cell &operator=(const Cell &src)¶
Assignment operator.
Called if one cell is assigned to another. Remember to change both assignment operator and copy constructor when adding new attributes to Cell.
-
inline int Colour(void) const¶
Returns the cell colour.
-
inline int SetColour(const int new_colour)¶
Set color of this cell to new_colour, irrespective of type.
-
inline void PrintInertia(void)¶
Debugging function used to print the cell’s current inertia tensor (as used for calculations of the length )
-
inline int Sigma() const¶
Returns the cell’s cell identity number.
-
inline int SetTargetArea(const int new_area)¶
Sets the target area of the cell.
-
inline void Apoptose()¶
Sends the current cell into apoptosis.
-
inline int IncrementTargetArea()¶
Decrement the cell’s target area by one unit.
-
inline int DecrementTargetArea()¶
Increment the cell’s target area by one unit.
-
inline int TimesDivided(void) const¶
Returns a counter keeping track of the number of divisions.
-
inline int DateOfBirth(void) const¶
Returns Monte Carlo Step (MCS) when this cell originated.
-
inline int ColourOfBirth(void) const¶
Returns the cell type at the time of birth.
-
inline double *SetGrad(double *g)¶
Set the current gradient of the cell to g. Currently not in use.
-
inline const double *GetGrad(void) const¶
Returns the cell’s measured gradient. Currently not in use.
-
inline double GradX() const¶
Returns the cell’s measured gradient. Currently not in use.
-
inline double GradY() const¶
Returns the cell’s measured gradient. Currently not in use.
-
inline double *AddToGrad(double *g)¶
Currently not in use (remove?)
-
inline void ClearGrad(void)¶
Currently not in use (remove?)
-
void MeasureCellSize(Cell &c)¶
After introducing a new Cell (e.g. with GrowInCell) call this function to set the moments and areas right.
-
inline int IncrementAdhesiveArea(int increment)¶
Increments the cell’s actual adhesive area by 1 unit.
-
inline int DecrementAdhesiveArea(int decrement)¶
Decrement the cell’s actual adhesive area by 1 unit.
Public Static Functions
-
static void ClearJ(void)¶
Clears the table of J’s.
This is only important for a feature called “DynamicJ’s”, where J-values depend on internal states of the cells (such as a genetic network; see e.g. Hogeweg et al. 2000). The current version of TST does not include such functionality.
-
static inline int MaxSigma()¶
Returns the maximum cell identity number in the Dish. This would normally be the number of cells in the Dish, although the number includes apoptosed cells.
-
static inline int SetJ(int t1, int t2, int val)¶
Sets bond energy J between cell type t1 and t2 to val.
-
inline Cell(const Dish &who, int settau = 1)¶
-
class PDE¶
Public Functions
-
PDE(const int layers, const int sizex, const int sizey)¶
Constructor for PDE object containing arbitrary number of planes.
PRIVATE
-
void Plot(Graphics *g, const int layer = 0)¶
Plots one layer of the PDE plane to a Graphics window.
- Parameters:
g – Graphics window.
layer – The PDE plane to be plotted. Default layer 0.
-
void Plot(Graphics *g, CellularPotts *cpm, const int layer = 0)¶
Plots one layer of the PDE to a Graphics window, but not over the cells.
- Parameters:
g – Graphics window.
cpm – CellularPotts object containing the cells.
layer – The PDE plane to be plotted. Default layer 0.
-
void ContourPlot(Graphics *g, int layer = 0, int colour = 1)¶
Plots the PDE field using contour lines.
- Parameters:
g – Graphics window.
layer – The PDE plane to be plotted. Default layer 0.
colour – Color to use for the contour lines, as defined in the “default.ctb” color map file, which should be in the same directory as the executable. Default color 1 (black in the default color map).
-
void SetSpeciesName(int l, const char *name)¶
Set the.
- Parameters:
name – of the species in layer
l –
-
inline PDEFIELD_TYPE get_PDEvars(const int layer, const int x, const int y) const¶
Returns the value of grid point x,y of PDE plane “layer”.
Warning, no range checking done.
- Parameters:
layer – the PDE plane to probe.
x, y – grid point to probe.
-
inline void setValue(const int layer, const int x, const int y, const PDEFIELD_TYPE value)¶
Sets grid point x,y of PDE plane “layer” to value “value”.
- Parameters:
layer – PDE plane.
x, y – grid point
value – new contents
-
inline void addtoValue(const int layer, const int x, const int y, const PDEFIELD_TYPE value)¶
Adds a number to a PDE grid point.
- Parameters:
layer – PDE plane.
x, y – grid point
value – value to add
-
inline PDEFIELD_TYPE Max(int l)¶
Gets the maximum value of PDE layer l.
- Parameters:
l – layer
- Returns:
Maximum value in layer l.
-
inline PDEFIELD_TYPE Min(int l)¶
Returns the minimum value of PDE layer l.
- Parameters:
l – layer
- Returns:
Minimum value in layer l.
-
void InitialiseAgeLayer(int l, double value, CellularPotts *cpm)¶
Carry out $n$ diffusion steps for all PDE planes. We use a forward Euler method here. Can be replaced for better algorithm. Function for the Act model. The whole field is initialised, usually with 0.
-
void NoFluxBoundaries(void)¶
Implementation of no-flux boundaries.
Called internally (optionally) by Diffuse().
-
void AbsorbingBoundaries(void)¶
Implementation of absorbing boundaries.
Called internally (optionally) by Diffuse().
-
void PeriodicBoundaries(void)¶
Implementation of periodic boundaries. Called internally (optionally) by Diffuse().
-
void InitialiseDiffusionCoefficients(CellularPotts *cpm)¶
Intialisation of diffusion coefficients.
- Parameters:
cpm – CellularPotts plane the PDE plane interacts with The initial diffusion coefficients may be space dependent on the cpm configuration
-
void InitialisePDE(CellularPotts *cpm)¶
Intialisation of PDE variables.
- Parameters:
cpm – CellularPotts plane the PDE plane interacts with Initial conditions conditions for the PDE should be given here.
-
void DerivativesPDE(CellularPotts *cpm, PDEFIELD_TYPE *derivs, int x, int y)¶
Derivatives of PDE variables.
- Parameters:
cpm – CellularPotts plane the PDE plane interacts with You should implement this member function as part of your main simulation code. See for an example vessel.cpp.
- Returns:
Derivatives at pixel (x,y)
-
void ForwardEulerStep(int repeat, CellularPotts *cpm)¶
Do a single forward Euler step to solve the ODE.
- Parameters:
repeat – Number of steps. We solve with a simple forward Euler solver. Ths can be replaced with alternative ODE solvers.
-
void Diffuse(int repeat)¶
Carry out $n$ diffusion steps for all PDE planes.
We use a forward Euler method here. Can be replaced for better algorithm.
Time step dt, space step dx, diffusion coefficient diff_coeff and boundary conditions (bool periodic_boundary) are set as global parameters in a parameter file using class Parameter.
- Parameters:
repeat – Number of steps.
-
void ReactionDiffusion(CellularPotts *cpm)¶
Do a single reaction diffusion step based on the given PDE derivatives.
-
void Secrete(CellularPotts *cpm)¶
Reaction and interaction of CPM plane with PDE planes.
- Parameters:
cpm – CellularPotts plane the PDE plane interacts with You should implement this member function as part of your main simulation code. See for an example vessel.cpp. This method is slightly faster than the general PDE solver.
-
inline double TheTime(void) const¶
Returns cumulative “simulated” time, i.e. number of time steps * dt.
-
double GetChemAmount(const int layer = -1)¶
Returns summed amount of chemical in PDE plane “layer”.
- Parameters:
layer – The PDE plane of which to sum the chemicals. layer=-1 (default) returns the summed amount of chemical in all planes.
-
void GradC(int layer = 0, int first_grad_layer = 1)¶
Calculates the first and second order gradients, i.e. gradx, grady, gradxx, gradxy and gradyy and puts them in the next three chemical fields. Not currently used and might need some redoing. Make sure you have allocated sufficient fields (this method generates five planes).
- Parameters:
layer – PDE plane of which to calculate the gradients (default 0)
first_grad_layer – first plane of five in which to write the results (default 1).
-
void PlotVectorField(Graphics &g, int stride, int linelength, int first_grad_layer = 1)¶
Plots a field of the first order gradients, i.e. gradx and grady; assumes you have called GradC before. Not currently used and might need some redoing.
- Parameters:
g – Graphics window
stride – Number of grid points between vectors (drawn as lines, currently.
linelength – Length of vector lines, in pixels.
first_grad_layer – first plane of two which contain the calculated gradients (default 1).
-
void InitLinearYGradient(int spec, double conc_top, double conc_bottom)¶
Initialise a linear gradient of the PDE variables in the Y direction.
- Parameters:
spec – The layer that will be initialised
conc_top – Concentration at the top of the matrix
conc_bottom – Concentration at the bottom of the matrix
-
bool plotPos(int x, int y, Graphics *graphics, int layer)¶
Plots the PDE variables of a pixel.
- Parameters:
x – x-coordinate
y – y-coordinate
graphics – Graphics interface
layer – Layer that will be displayed
-
void InitialiseCuda()¶
allocate memory required for the CUDA reaction-diffusion solver To use this CUDA solver, the following steps must be taken. -Make sure you have an Nvidia GPU -Install CUDA 12.x or higher (contains required cuSparse version for cusparseDgtsvInterleavedBatch and cusparseSgtsvInterleavedBatch) -Implement your derivatives function in device void DerivativesPDE in pde.cu -Enable CUDA by changing the USECUDA flag in Tissue_Similation_Toolkit.pri -Recompile your code base by performing ‘make clean’ and ‘qmake’ -Specify the desired number of cores and threads per core in the parameter file -Enable CUDA by using ‘usecuda = true’ in the parameter file
-
void InitialisePDEvars(CellularPotts *cpm, int *celltypes)¶
Intialise the PDE variables The variables will be used to solve the reaction diffusion equation.
- Parameters:
cpm – Initialisation may depend on the CPM configuration
celltypes – Initalisation may depend on the celltypes This initialisation is used for both the CPU and CUDA solver
-
void cuPDEsteps(CellularPotts *cpm, int repeats)¶
Do a single reaction diffusion step on CUDA of size dt A reaction-diffusion step of time dt is performed on CUDA. The reaction part is solved with the Forward Euler method on CUDA. The diffusion is solved with the alternating directions implicit (ADI) method. Communication between host and device is only necessary before and after the reaction-diffusion step. The PDEvars vector is updated accordingly. PDEvars need to be initialised with InitialisePDEvars. The derivatives must be described in device void DerivativesPDE in pde.cu, rather than the PDE derivatives function from pde.hpp A single reaction-diffusion step is structured as follows:
Perform an Forward Euler steps of size dt/2 in increments of ddt
Perform a horizontal ADI sweep of size dt/2
Perform an Forward Euler steps of size dt/2 in increments of ddt
Perform a horizontal ADI sweep of size dt/2
- Parameters:
cpm – The current CPM configuration
repeats – Number of reaction-diffusion steps that are performed
-
void cuODEstep(void)¶
Perform a single ODE step Currently this is done with a forward Euler solver, but other solvers may be implemented in a straight forward way.
-
void cuHorizontalADIstep(void)¶
Perform the horizontal alternating directions implicit method step This uses an interleaved format for solving which is taken care of by the initialisation functions.
-
void cuVerticalADIstep(void)¶
Perform the vertical alternating directions implicit method step This uses an interleaved format for solving which is taken care of by the initialisation functions.
-
PDE(const int layers, const int sizex, const int sizey)¶
ECM and adhesion code¶
-
class AdhesionIndex¶
Tracks location of adhesions in the ECM grid.
This class provides the adhesion particles and their bonds and angle constraints per ECM pixel, in a format that allows for efficient force calculations. Bonds or angle constraints involving other particles that are themselves adhesions or that are excluded are ignored.
Public Functions
-
void rebuild(ECMBoundaryState const &ecm_boundary)¶
Rebuild the cached data to match the ECM boundary again.
This must be called after every change to the ECM, except for changes that only move adhesions from one pixel to another, for which move_adhesions() should be used because it’s more efficient.
- Parameters:
ecm_boundary – The current state of the ECM boundary
-
std::vector<AdhesionWithEnvironment> const &get_adhesions(PixelPos pixel) const¶
Get adhesions at a given pixel.
Note that this function returns a reference to a vector. This reference will be invalidated by any subsequent call to update() or move_adhesions() on this object.
- Parameters:
pixel – Pixel for which to get adhesions.
-
void move_adhesions(PixelPos from, PixelPos to)¶
Move adhesions from one pixel to another.
This modifies only the cache, the ECM needs to be updated separately.
- Parameters:
from – Pixel to move adhesions from
to – Pixel to move adhesions to
-
void remove_adhesions(PixelPos pixel)¶
Remove adhesions at the given pixel.
This modifies only the cache, the ECM needs to be updated separately.
- Parameters:
pixel – Pixel to move adhesions from
-
CellECMInteractions get_cell_ecm_interactions() const¶
Get accumulated changes to the adhesions.
AdhesionIndex keeps track of any changes to the adhesions it applies. This function returns the accumulated changes.
-
void reset_cell_ecm_interactions()¶
Reset the adhesion change administration.
This clears the recorded adhesion change history.
-
void rebuild(ECMBoundaryState const &ecm_boundary)¶
-
class AdhesionMover¶
Decide how adhesions will move during a copy.
This class implements the logic required to evaluate a copy attempt that requires moving one or more adhesion particles along with or out of the way of the copied pixel.
Public Functions
-
AdhesionMover(CellularPotts const &ca)¶
Construct an AdhesionMover object.
- Parameters:
ca – The CPM to work with.
-
double move_dh(PixelPos from, PixelPos to, AdhesionDisplacements &displacements) const¶
Compute the amount of work involved in an attempted copy.
A Cellular Potts copy attempt from and/or to a grid cell that has adhesions may move those adhesions along, dragging along the ECM in turn. This involves work (a change in the overall energy of the system), the amount of which is calculated by this function for a given potential move.
Note that there are multiple options for whether and where to move the adhesions on a copy attempt, in particular for retractions. This function returns the work required along with the corresponding chosen displacement. If the copy attempt succeeds, then the latter should be passed back to commit_move() to ensure the correct displacement is made.
- Parameters:
from – Pixel to be copied from
to – Pixel to copied to
displacements – Out parameter returning the chosen adhesion displacements
- Returns:
The work required.
-
void commit_move(PixelPos source_pixel, PixelPos target_pixel, AdhesionDisplacements const &displacements)¶
Update the adhesions following a move.
This modifies the associated ECM as well as the cache, so update() does not need to be called afterwards.
This will not do a simultaneous swap if the displacements are in opposite directions aligned with the move. That should never happen, as it requires attempting a copy within a cell which would have been short-circuited long before this got called.
- Parameters:
from – Pixel that was copied from
to – Pixel that was copied to
displacements – Adhesion displacements as previously returned by move_dh()
-
CellECMInteractions get_cell_ecm_interactions() const¶
Get accumulated changes to the adhesions.
AdhesionMover keeps track of any changes to the adhesions it applies. This function returns the accumulated changes.
-
void reset_cell_ecm_interactions()¶
Reset the adhesion change administration.
This clears the recorded adhesion change history.
-
void update(ECMBoundaryState const &ecm_boundary)¶
Update the internal administration after an update to the ECM.
AdhesionMover keeps a partial copy of the ECM internally for faster access. It keeps this up to date when any adhesions are moved via commit_move, but if external changes are made to the ECM then this function must be called to get the AdhesionMover back in sync.
- Parameters:
ecm_boundary – ECM boundary state to update from
-
AdhesionMover(CellularPotts const &ca)¶
-
class AdhesionDisplacements¶
Describes where the adhesions go during a copy attempt.
Public Functions
-
AdhesionDisplacements()¶
Create an object where no adhesions move.
-
AdhesionDisplacements(PixelDisplacement source, PixelDisplacement target)¶
Create an object representing the given displacements.
- Parameters:
source – Displacement of adhesions in the source pixel.
target – Displacement of adhesions in the target pixel.
Public Members
-
PixelDisplacement source¶
Location relative to source pixel where adhesions go.
-
PixelDisplacement target¶
Location relative to target pixel where adhesions go.
Public Static Attributes
-
static const PixelDisplacement annihilated¶
Special value signalling that the adhesions got annihilated.
-
AdhesionDisplacements()¶