-
Notifications
You must be signed in to change notification settings - Fork 4
Expand file tree
/
Copy pathitkBioCellularAggregate.h
More file actions
197 lines (155 loc) · 5.8 KB
/
itkBioCellularAggregate.h
File metadata and controls
197 lines (155 loc) · 5.8 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#ifndef itkBioCellularAggregate_h
#define itkBioCellularAggregate_h
#include "itkDefaultDynamicMeshTraits.h"
#include "itkMesh.h"
#include "itkImage.h"
#include "itkBioCell.h"
#include "itkPolygonCell.h"
#include <iostream>
#include <vector>
namespace itk
{
namespace bio
{
/** \class CellularAggregate
* \brief Base class for different types of cellular groups,
* including bacterial colonies and pluricellular organisms
*
* This class represents an aggregation of bio::Cell objects.
*
* \ingroup ITKBioCell
*/
template <unsigned int NSpaceDimension = 3>
class ITK_TEMPLATE_EXPORT CellularAggregate : public CellularAggregateBase
{
public:
ITK_DISALLOW_COPY_AND_MOVE(CellularAggregate);
/** Standard class type alias. */
using Self = CellularAggregate;
using Superclass = CellularAggregateBase;
using Pointer = SmartPointer<Self>;
using ConstPointer = SmartPointer<const Self>;
/*** Run-time type information (and related methods). */
itkTypeMacro(CellularAggregate, CellularAggregateBase);
/** Method for creation through the object factory. */
itkNewMacro(Self);
static constexpr unsigned int SpaceDimension = NSpaceDimension;
/*** Type to be used for data associated with each point in the mesh. */
using BioCellType = Cell<NSpaceDimension>;
using PointPixelType = BioCellType *;
using CellPixelType = double;
/** Mesh Traits */
using MeshTraits = DefaultDynamicMeshTraits<PointPixelType, // PixelType
NSpaceDimension, // Points Dimension
NSpaceDimension, // Max.Topological Dimension
double, // Type for coordinates
double, // Type for interpolation
CellPixelType // Type for values in the cells
>;
/** Mesh Traits */
using MeshType = Mesh<PointPixelType, NSpaceDimension, MeshTraits>;
/** Mesh Associated types */
using MeshPointer = typename MeshType::Pointer;
using MeshConstPointer = typename MeshType::ConstPointer;
using PointType = typename MeshType::PointType;
using VectorType = typename BioCellType::VectorType;
using PointsContainer = typename MeshType::PointsContainer;
using PointDataContainer = typename MeshType::PointDataContainer;
using VoronoiRegionsContainer = typename MeshType::CellsContainer;
using PointsIterator = typename PointsContainer::Iterator;
using CellsIterator = typename PointDataContainer::Iterator;
using VoronoiIterator = typename VoronoiRegionsContainer::Iterator;
using PointsConstIterator = typename PointsContainer::ConstIterator;
using CellsConstIterator = typename PointDataContainer::ConstIterator;
using VoronoiConstIterator = typename VoronoiRegionsContainer::ConstIterator;
using CellAutoPointer = typename MeshType::CellAutoPointer;
/** Voronoi region around a bio::Cell */
using CellInterfaceType = CellInterface<typename MeshType::CellPixelType, typename MeshType::CellTraits>;
using VoronoiRegionType = PolygonCell<CellInterfaceType>;
using VoronoiRegionAutoPointer = typename VoronoiRegionType::SelfAutoPointer;
/** Convenient type alias. */
using ImagePixelType = float;
using SubstrateType = Image<ImagePixelType, NSpaceDimension>;
using SubstratePointer = typename SubstrateType::Pointer;
using SubstrateValueType = ImagePixelType;
using SubstratesVector = std::vector<SubstratePointer>;
public:
unsigned int
GetNumberOfCells() const;
static unsigned int
GetDimension()
{
return SpaceDimension;
}
void
SetGrowthRadiusLimit(double value);
void
SetGrowthRadiusIncrement(double value);
itkGetModifiableObjectMacro(Mesh, MeshType);
virtual void
AdvanceTimeStep();
virtual void
SetEgg(BioCellType * cell, const PointType & position);
virtual void
Add(CellBase * cell);
virtual void
Add(CellBase * cell, const VectorType & perturbation);
void
Add(CellBase * cellA, CellBase * cellB, double perturbationLength) override;
void
Remove(CellBase * cell) override;
virtual void
GetVoronoi(IdentifierType cellId, VoronoiRegionAutoPointer &) const;
void
DumpContent(std::ostream & os) const;
virtual void
AddSubstrate(SubstrateType * substrate);
virtual SubstratesVector &
GetSubstrates();
SubstrateValueType
GetSubstrateValue(IdentifierType cellId, unsigned int substrateId) const override;
virtual void
KillAll();
protected:
CellularAggregate();
~CellularAggregate() override;
void
PrintSelf(std::ostream & os, Indent indent) const override;
virtual void
ComputeForces();
virtual void
UpdatePositions();
virtual void
ComputeClosestPoints();
virtual void
ClearForces();
private:
MeshPointer m_Mesh;
SubstratesVector m_Substrates;
double m_FrictionForce;
SizeValueType m_Iteration;
SizeValueType m_ClosestPointComputationInterval;
};
} // end namespace bio
} // end namespace itk
#ifndef ITK_MANUAL_INSTANTIATION
# include "itkBioCellularAggregate.hxx"
#endif
#endif