Shapeworks Studio  2.1
Shape analysis software suite
itkPSMShapeEntropyFunctionTest.cxx
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 #include <iostream>
19 #include "itkPSMShapeEntropyFunction.h"
20 #include "itkPSMParticleSystem.h"
21 #include "itkPSMRegionNeighborhood.h"
22 #include "itkPSMRegionDomain.h"
23 #include "itkMacro.h"
24 
26 int itkPSMShapeEntropyFunctionTest(int, char* [] )
27 {
28  bool passed = true;
29  std::string errstring = "";
30 
31  try
32  {
33  typedef itk::Point<double, 3> PointType;
34  itk::PSMShapeEntropyFunction<3>::Pointer p = itk::PSMShapeEntropyFunction<3>::New();
35  itk::PSMShapeMatrixAttribute<double, 3>::Pointer sm = itk::PSMShapeMatrixAttribute<double, 3>::New();
36  sm->SetDomainsPerShape(1);
37  p->SetShapeMatrix(sm);
38 
39  itk::PSMParticleSystem<3>::Pointer ps = itk::PSMParticleSystem<3>::New();
40  ps->RegisterAttribute(sm);
41 
42  // fill in particle system with points
43  itk::PSMRegionNeighborhood<3>::Pointer nb1 = itk::PSMRegionNeighborhood<3>::New();
44  itk::PSMRegionNeighborhood<3>::Pointer nb2 = itk::PSMRegionNeighborhood<3>::New();
45  itk::PSMRegionDomain<3>::Pointer d1 = itk::PSMRegionDomain<3>::New();
46  itk::PSMRegionDomain<3>::Pointer d2 = itk::PSMRegionDomain<3>::New();
47 
48  // Define bounding box
49  const unsigned int SZ = 100;
50  PointType ptl, ptu;
51  ptl[0] = 0.0f; ptl[1] = 0.0f; ptl[2] = 0.0f;
52  ptu[0] = static_cast<double>(SZ); ptu[1] = static_cast<double>(SZ); ptu[2] = static_cast<double>(SZ);
53 
54  // Add domains and neighborhoods
55  d1->SetRegion(ptl, ptu);
56  d2->SetRegion(ptl, ptu);
57  ps->AddDomain(d1);
58  ps->AddDomain(d2);
59  ps->SetNeighborhood(0, nb1);
60  ps->SetNeighborhood(1, nb2);
61 
62  // Add position
63  PointType pt;
64  for (unsigned int i = 0; i < SZ; i++)
65  {
66  pt[0] = static_cast<double>(i) + 0.1f;
67  pt[1] = static_cast<double>(i) + 0.2f;
68  pt[2] = static_cast<double>(i) + 0.3f;
69  ps->AddPosition(pt, 0);
70  ps->AddPosition(pt, 1);
71  }
72 
73  // Evaluate
74  p->BeforeIteration();
75  double energy00 = p->Energy(0, 0, ps);
76  double energy10 = p->Energy(1, 0, ps);
77  double energy01 = p->Energy(0, 1, ps);
78  double energy11 = p->Energy(1, 1, ps);
79  if (energy00 != energy01 || energy10 != energy11)
80  {
81  passed = false;
82  errstring += std::string("Energy method failed. ");
83  }
84  }
85  catch(itk::ExceptionObject &e)
86  {
87  errstring = "ITK exception with description: " + std::string(e.GetDescription())
88  + std::string("\n at location:") + std::string(e.GetLocation())
89  + std::string("\n in file:") + std::string(e.GetFile());
90  passed = false;
91  }
92  catch(...)
93  {
94  errstring = "Unknown exception thrown";
95  passed = false;
96  }
97 
98  if (passed)
99  {
100  std::cout << "All tests passed" << std::endl;
101  return EXIT_SUCCESS;
102  }
103  else
104  {
105  std::cout << "Test failed with the following error:" << std::endl;
106  std::cout << errstring << std::endl;
107  return EXIT_FAILURE;
108  }
109 }
Each column describes a shape. A shape may be composed of m_DomainsPerShape domains (default 1)...
void AddDomain(DomainType *, int threadId=0)
void SetNeighborhood(unsigned int, NeighborhoodType *, int threadId=0)
void RegisterAttribute(PSMAttribute< VDimension > *)
const PointType & AddPosition(const PointType &, unsigned int d=0, int threadId=0)
void SetShapeMatrix(ShapeMatrixType *s)
A facade class that manages interactions with a particle system.
void SetRegion(const PointType &l, const PointType &u)