PRISMS-PF  v2.1
test_SimplifiedGrainRepresentation.h
Go to the documentation of this file.
1 #include "../../include/SimplifiedGrainRepresentation.h"
2 
3 template <int dim,typename T>
5 
6  char buffer[100];
7 
8  std::cout << "\nTesting 'SimplifiedGrainRepresentation'... " << std::endl;
9 
10  bool pass = true;
11 
12  // Subtest 1 (single element)
13  {
14  GrainSet<dim> test_grain_set;
15 
16  test_grain_set.setGrainIndex(2);
17 
18  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
19  {dealii::Point<dim> p(0.0,0.75); vertex_set[0] = p;}
20  {dealii::Point<dim> p(0.25,0.75); vertex_set[1] = p;}
21  {dealii::Point<dim> p(0.0,1.0); vertex_set[2] = p;}
22  {dealii::Point<dim> p(0.25,1.0); vertex_set[3] = p;}
23 
24  test_grain_set.addVertexList(vertex_set);
25 
26  SimplifiedGrainRepresentation<dim> simplified_grain_representation(test_grain_set);
27 
28  std::cout << "Centroid: " << simplified_grain_representation.getCenter() << std::endl;
29  std::cout << "Radius: " << simplified_grain_representation.getRadius() << std::endl;
30 
31  bool result = false;
32  if ( (std::abs(simplified_grain_representation.getCenter()(0) - 0.125) < 1.0e-10) and (std::abs(simplified_grain_representation.getCenter()(1) - 0.875) < 1.0e-10) and (std::abs(simplified_grain_representation.getRadius() - 0.125*std::sqrt(2.0)) < 1.0e-10)){
33  result = true;
34  }
35 
36  sprintf (buffer, "Subtest 1 result for 'SimplifiedGrainRepresentation': %u\n", result);
37  std::cout << buffer;
38  pass = pass and result;
39  }
40 
41  // Subtest 2 (L-shaped grain)
42  {
43  GrainSet<dim> test_grain_set;
44 
45  test_grain_set.setGrainIndex(2);
46 
47  {
48  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
49  {dealii::Point<dim> p(0.0,0.75); vertex_set[0] = p;}
50  {dealii::Point<dim> p(0.25,0.75); vertex_set[1] = p;}
51  {dealii::Point<dim> p(0.0,1.0); vertex_set[2] = p;}
52  {dealii::Point<dim> p(0.25,1.0); vertex_set[3] = p;}
53 
54  test_grain_set.addVertexList(vertex_set);
55  }
56  {
57  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
58  {dealii::Point<dim> p(0.0,0.5); vertex_set[0] = p;}
59  {dealii::Point<dim> p(0.25,0.5); vertex_set[1] = p;}
60  {dealii::Point<dim> p(0.0,0.75); vertex_set[2] = p;}
61  {dealii::Point<dim> p(0.25,0.75); vertex_set[3] = p;}
62 
63  test_grain_set.addVertexList(vertex_set);
64  }
65  {
66  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
67  {dealii::Point<dim> p(0.25,0.75); vertex_set[0] = p;}
68  {dealii::Point<dim> p(0.5,0.75); vertex_set[1] = p;}
69  {dealii::Point<dim> p(0.25,1.0); vertex_set[2] = p;}
70  {dealii::Point<dim> p(0.5,1.0); vertex_set[3] = p;}
71 
72  test_grain_set.addVertexList(vertex_set);
73  }
74 
75  SimplifiedGrainRepresentation<dim> simplified_grain_representation(test_grain_set);
76 
77  std::cout << "Centroid: " << simplified_grain_representation.getCenter() << std::endl;
78  std::cout << "Radius: " << simplified_grain_representation.getRadius() << std::endl;
79 
80  bool result = false;
81 
82  double centroid_x = (0.125*3.0+0.25)/3.0;
83  double centroid_y = (0.875*3.0-0.25)/3.0;
84  double radius = std::sqrt(dealii::Utilities::fixed_power<2>( centroid_x - 0.5) + dealii::Utilities::fixed_power<2>( centroid_y - 1.0));
85 
86  if ( (std::abs(simplified_grain_representation.getCenter()(0) - centroid_x) < 1.0e-10) and (std::abs(simplified_grain_representation.getCenter()(1) - centroid_y) < 1.0e-10) and (std::abs(simplified_grain_representation.getRadius() - radius) < 1.0e-10)){
87  result = true;
88  }
89 
90  sprintf (buffer, "Subtest 2 result for 'SimplifiedGrainRepresentation': %u\n", result);
91  std::cout << buffer;
92  pass = pass and result;
93  }
94 
95  sprintf (buffer, "Test result for 'SimplifiedGrainRepresentation': %u\n", pass);
96  std::cout << buffer;
97 
98  return pass;
99 }
dealii::Point< dim > getCenter() const
void addVertexList(std::vector< dealii::Point< dim >> _vertices)
Definition: FloodFiller.h:48
void setGrainIndex(unsigned int _grain_index)
Definition: FloodFiller.h:28