ITK
4.1.0
Insight Segmentation and Registration Toolkit
|
00001 /*========================================================================= 00002 * 00003 * Copyright Insight Software Consortium 00004 * 00005 * Licensed under the Apache License, Version 2.0 (the "License"); 00006 * you may not use this file except in compliance with the License. 00007 * You may obtain a copy of the License at 00008 * 00009 * http://www.apache.org/licenses/LICENSE-2.0.txt 00010 * 00011 * Unless required by applicable law or agreed to in writing, software 00012 * distributed under the License is distributed on an "AS IS" BASIS, 00013 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. 00014 * See the License for the specific language governing permissions and 00015 * limitations under the License. 00016 * 00017 *=========================================================================*/ 00018 #ifndef __itkBioCellularAggregate_h 00019 #define __itkBioCellularAggregate_h 00020 00021 #include "itkDefaultDynamicMeshTraits.h" 00022 #include "itkMesh.h" 00023 #include "itkImage.h" 00024 #include "itkBioCell.h" 00025 #include "itkPolygonCell.h" 00026 00027 #include <iostream> 00028 #include <vector> 00029 00030 namespace itk 00031 { 00032 namespace bio 00033 { 00042 template< unsigned int NSpaceDimension = 3 > 00043 class ITK_EXPORT CellularAggregate:public CellularAggregateBase 00044 { 00045 public: 00047 typedef CellularAggregate Self; 00048 typedef CellularAggregateBase Superclass; 00049 typedef SmartPointer< Self > Pointer; 00050 typedef SmartPointer< const Self > ConstPointer; 00051 00052 /*** Run-time type information (and related methods). */ 00053 itkTypeMacro(BioCellularAggregate, CellularAggregateBase); 00054 00056 itkNewMacro(Self); 00057 00058 itkStaticConstMacro(SpaceDimension, unsigned int, NSpaceDimension); 00059 00060 /*** Type to be used for data associated with each point in the mesh. */ 00061 typedef Cell< NSpaceDimension > BioCellType; 00062 typedef BioCellType * PointPixelType; 00063 typedef double CellPixelType; 00064 00066 typedef DefaultDynamicMeshTraits< 00067 PointPixelType, // PixelType 00068 NSpaceDimension, // Points Dimension 00069 NSpaceDimension, // Max.Topological Dimension 00070 double, // Type for coordinates 00071 double, // Type for interpolation 00072 CellPixelType // Type for values in the cells 00073 > MeshTraits; 00074 00076 typedef Mesh< PointPixelType, 00077 NSpaceDimension, 00078 MeshTraits > MeshType; 00079 00081 typedef typename MeshType::Pointer MeshPointer; 00082 typedef typename MeshType::ConstPointer MeshConstPointer; 00083 typedef typename MeshType::PointType PointType; 00084 typedef typename BioCellType::VectorType VectorType; 00085 00086 typedef typename MeshType::PointsContainer PointsContainer; 00087 typedef typename MeshType::PointDataContainer PointDataContainer; 00088 typedef typename MeshType::CellsContainer VoronoiRegionsContainer; 00089 typedef typename PointsContainer::Iterator PointsIterator; 00090 typedef typename PointDataContainer::Iterator CellsIterator; 00091 typedef typename VoronoiRegionsContainer::Iterator VoronoiIterator; 00092 typedef typename PointsContainer::ConstIterator PointsConstIterator; 00093 typedef typename PointDataContainer::ConstIterator CellsConstIterator; 00094 typedef typename VoronoiRegionsContainer::ConstIterator VoronoiConstIterator; 00095 typedef typename MeshType::CellAutoPointer CellAutoPointer; 00096 00098 typedef CellInterface< 00099 typename MeshType::CellPixelType, 00100 typename MeshType::CellTraits > CellInterfaceType; 00101 00102 typedef PolygonCell< CellInterfaceType > VoronoiRegionType; 00103 typedef typename VoronoiRegionType::SelfAutoPointer VoronoiRegionAutoPointer; 00104 00106 typedef float ImagePixelType; 00107 typedef Image< ImagePixelType, NSpaceDimension > SubstrateType; 00108 typedef typename SubstrateType::Pointer SubstratePointer; 00109 typedef ImagePixelType SubstrateValueType; 00110 typedef std::vector< SubstratePointer > SubstratesVector; 00111 public: 00112 unsigned int GetNumberOfCells(void) const; 00113 00114 static unsigned int GetDimension() { return SpaceDimension; } 00115 00116 void SetGrowthRadiusLimit(double value); 00117 00118 void SetGrowthRadiusIncrement(double value); 00119 00120 itkGetObjectMacro(Mesh, MeshType); 00121 itkGetConstObjectMacro(Mesh, MeshType); 00122 00123 virtual void AdvanceTimeStep(void); 00124 00125 virtual void SetEgg(BioCellType *cell, const PointType & position); 00126 00127 virtual void Add(CellBase *cell); 00128 00129 virtual void Add(CellBase *cell, const VectorType & perturbation); 00130 00131 virtual void Add(CellBase *cellA, CellBase *cellB, double perturbationLength); 00132 00133 virtual void Remove(CellBase *cell); 00134 00135 virtual void GetVoronoi(IdentifierType cellId, VoronoiRegionAutoPointer &) const; 00136 00137 void DumpContent(std::ostream & os) const; 00138 00139 virtual void AddSubstrate(SubstrateType *substrate); 00140 00141 virtual SubstratesVector & GetSubstrates(void); 00142 00143 virtual SubstrateValueType GetSubstrateValue(IdentifierType cellId, 00144 unsigned int substrateId) const; 00145 00146 virtual void KillAll(void); 00147 00148 protected: 00149 CellularAggregate(); 00150 virtual ~CellularAggregate(); 00151 CellularAggregate(const Self &); 00152 void operator=(const Self &); 00153 00154 void PrintSelf(std::ostream & os, Indent indent) const; 00155 00156 virtual void ComputeForces(void); 00157 00158 virtual void UpdatePositions(void); 00159 00160 virtual void ComputeClosestPoints(void); 00161 00162 virtual void ClearForces(void); 00163 00164 private: 00165 00166 MeshPointer m_Mesh; 00167 SubstratesVector m_Substrates; 00168 double m_FrictionForce; 00169 SizeValueType m_Iteration; 00170 SizeValueType m_ClosestPointComputationInterval; 00171 }; 00172 } // end namespace bio 00173 } // end namespace itk 00174 00175 #ifndef ITK_MANUAL_INSTANTIATION 00176 #include "itkBioCellularAggregate.hxx" 00177 #endif 00178 00179 #endif 00180