ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkBioCellularAggregate.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkBioCellularAggregate_h
19 #define itkBioCellularAggregate_h
20 
22 #include "itkMesh.h"
23 #include "itkImage.h"
24 #include "itkBioCell.h"
25 #include "itkPolygonCell.h"
26 
27 #include <iostream>
28 #include <vector>
29 
30 namespace itk
31 {
32 namespace bio
33 {
42 template< unsigned int NSpaceDimension = 3 >
44 {
45 public:
51 
52 /*** Run-time type information (and related methods). */
53  itkTypeMacro(BioCellularAggregate, CellularAggregateBase);
54 
56  itkNewMacro(Self);
57 
58  itkStaticConstMacro(SpaceDimension, unsigned int, NSpaceDimension);
59 
60 /*** Type to be used for data associated with each point in the mesh. */
63  typedef double CellPixelType;
64 
67  PointPixelType, // PixelType
68  NSpaceDimension, // Points Dimension
69  NSpaceDimension, // Max.Topological Dimension
70  double, // Type for coordinates
71  double, // Type for interpolation
72  CellPixelType // Type for values in the cells
74 
76  typedef Mesh< PointPixelType,
77  NSpaceDimension,
79 
81  typedef typename MeshType::Pointer MeshPointer;
83  typedef typename MeshType::PointType PointType;
85 
89  typedef typename PointsContainer::Iterator PointsIterator;
90  typedef typename PointDataContainer::Iterator CellsIterator;
91  typedef typename VoronoiRegionsContainer::Iterator VoronoiIterator;
92  typedef typename PointsContainer::ConstIterator PointsConstIterator;
93  typedef typename PointDataContainer::ConstIterator CellsConstIterator;
94  typedef typename VoronoiRegionsContainer::ConstIterator VoronoiConstIterator;
96 
98  typedef CellInterface<
99  typename MeshType::CellPixelType,
101 
103  typedef typename VoronoiRegionType::SelfAutoPointer VoronoiRegionAutoPointer;
104 
106  typedef float ImagePixelType;
110  typedef std::vector< SubstratePointer > SubstratesVector;
111 
112 public:
113  unsigned int GetNumberOfCells() const;
114 
115  static unsigned int GetDimension() { return SpaceDimension; }
116 
117  void SetGrowthRadiusLimit(double value);
118 
119  void SetGrowthRadiusIncrement(double value);
120 
121  itkGetModifiableObjectMacro(Mesh, MeshType);
122 
123  virtual void AdvanceTimeStep();
124 
125  virtual void SetEgg(BioCellType *cell, const PointType & position);
126 
127  virtual void Add(CellBase *cell);
128 
129  virtual void Add(CellBase *cell, const VectorType & perturbation);
130 
131  virtual void Add(CellBase *cellA, CellBase *cellB, double perturbationLength) ITK_OVERRIDE;
132 
133  virtual void Remove(CellBase *cell) ITK_OVERRIDE;
134 
135  virtual void GetVoronoi(IdentifierType cellId, VoronoiRegionAutoPointer &) const;
136 
137  void DumpContent(std::ostream & os) const;
138 
139  virtual void AddSubstrate(SubstrateType *substrate);
140 
141  virtual SubstratesVector & GetSubstrates();
142 
144  unsigned int substrateId) const ITK_OVERRIDE;
145 
146  virtual void KillAll();
147 
148 protected:
150  virtual ~CellularAggregate();
151  CellularAggregate(const Self &);
152  void operator=(const Self &);
153 
154  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
155 
156  virtual void ComputeForces();
157 
158  virtual void UpdatePositions();
159 
160  virtual void ComputeClosestPoints();
161 
162  virtual void ClearForces();
163 
164 private:
165 
171 };
172 } // end namespace bio
173 } // end namespace itk
174 
175 #ifndef ITK_MANUAL_INSTANTIATION
176 #include "itkBioCellularAggregate.hxx"
177 #endif
178 
179 #endif
void SetGrowthRadiusIncrement(double value)
unsigned int GetNumberOfCells() const
Base class for different types of cellular groups, including bacterial colonies and pluricellular org...
virtual void ComputeForces()
Light weight base class for most itk classes.
MeshTraits::PointDataContainer PointDataContainer
Definition: itkMesh.h:157
virtual void Add(CellBase *cell)
Image< ImagePixelType, NSpaceDimension > SubstrateType
virtual void AdvanceTimeStep()
void SetGrowthRadiusLimit(double value)
CellInterface< typename MeshType::CellPixelType, typename MeshType::CellTraits > CellInterfaceType
MeshType::ConstPointer MeshConstPointer
virtual void GetVoronoi(IdentifierType cellId, VoronoiRegionAutoPointer &) const
DefaultDynamicMeshTraits< PointPixelType, NSpaceDimension, NSpaceDimension, double, double, CellPixelType > MeshTraits
MeshTraits::PointsContainer PointsContainer
Definition: itkMesh.h:152
virtual SubstratesVector & GetSubstrates()
MeshTraits::PointType PointType
Definition: itkMesh.h:151
Implements the N-dimensional mesh structure.
Definition: itkMesh.h:108
virtual void PrintSelf(std::ostream &os, Indent indent) const override
An abstract interface for cells.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
This class implements the minimal behavior of a biological cell.
Definition: itkBioCell.h:39
BioCellType::VectorType VectorType
SizeValueType IdentifierType
Definition: itkIntTypes.h:147
PolygonCell< CellInterfaceType > VoronoiRegionType
virtual void ComputeClosestPoints()
MeshType::PointDataContainer PointDataContainer
Base class for the CellularAggregates.
PointDataContainer::ConstIterator CellsConstIterator
PointDataContainer::Iterator CellsIterator
PointsContainer::ConstIterator PointsConstIterator
VoronoiRegionsContainer::Iterator VoronoiIterator
SmartPointer< const Self > ConstPointer
CellType::CellAutoPointer CellAutoPointer
Definition: itkMesh.h:190
PointsContainer::Iterator PointsIterator
virtual SubstrateValueType GetSubstrateValue(IdentifierType cellId, unsigned int substrateId) const override
Represents a polygon in a Mesh.
void DumpContent(std::ostream &os) const
virtual void UpdatePositions()
virtual void SetEgg(BioCellType *cell, const PointType &position)
Mesh< PointPixelType, NSpaceDimension, MeshTraits > MeshType
virtual void Remove(CellBase *cell) override
MeshTraits::CellsContainer CellsContainer
Definition: itkMesh.h:154
void operator=(const Self &)
Non-templated Base class from which the templated Cell classes will be derived.
VoronoiRegionType::SelfAutoPointer VoronoiRegionAutoPointer
MeshTraits::CellPixelType CellPixelType
Definition: itkMesh.h:128
Control indentation during Print() invocation.
Definition: itkIndent.h:49
SubstrateType::Pointer SubstratePointer
A simple structure that holds type information for a mesh and its cells.
MeshTraits::CellTraits CellTraits
Definition: itkMesh.h:153
MeshType::PointsContainer PointsContainer
VoronoiRegionsContainer::ConstIterator VoronoiConstIterator
std::vector< SubstratePointer > SubstratesVector
Templated n-dimensional image class.
Definition: itkImage.h:75
MeshType::CellsContainer VoronoiRegionsContainer
Cell< NSpaceDimension > BioCellType
MeshType::CellAutoPointer CellAutoPointer
static const unsigned int SpaceDimension
virtual void AddSubstrate(SubstrateType *substrate)