PRISMS-PF  v2.1
test_SimplifiedGrainManipulator.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 'SimplifiedGrainManipulator::transferGrainIds'... " << std::endl;
9 
10  bool pass = true;
11 
12  // Subtest 1 (two grains both needing reassignment)
13  {
14  // Create grain 0_old with center (0.5, 0.5), radius 0.5*sqrt(2), and order parameter 0
15  GrainSet<dim> test_grain_set_0;
16  {
17  test_grain_set_0.setGrainIndex(0);
18  test_grain_set_0.setOrderParameterIndex(0);
19 
20  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
21  {dealii::Point<dim> p(0.0,0.0); vertex_set[0] = p;}
22  {dealii::Point<dim> p(1.0,0.0); vertex_set[1] = p;}
23  {dealii::Point<dim> p(0.0,1.0); vertex_set[2] = p;}
24  {dealii::Point<dim> p(1.0,1.0); vertex_set[3] = p;}
25 
26  test_grain_set_0.addVertexList(vertex_set);
27  }
28 
29  SimplifiedGrainRepresentation<dim> simplified_grain_representation_0(test_grain_set_0);
30 
31  // Create grain 1_old with center (4.5, 0.5), radius 0.5*sqrt(2), and order parameter 1
32  GrainSet<dim> test_grain_set_1;
33  {
34  test_grain_set_1.setGrainIndex(1);
35  test_grain_set_1.setOrderParameterIndex(1);
36 
37  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
38  {dealii::Point<dim> p(4.0,0.0); vertex_set[0] = p;}
39  {dealii::Point<dim> p(5.0,0.0); vertex_set[1] = p;}
40  {dealii::Point<dim> p(4.0,1.0); vertex_set[2] = p;}
41  {dealii::Point<dim> p(5.0,1.0); vertex_set[3] = p;}
42 
43  test_grain_set_1.addVertexList(vertex_set);
44  }
45 
46  SimplifiedGrainRepresentation<dim> simplified_grain_representation_1(test_grain_set_1);
47 
48  // Create grain 0_new with center (0.6, 0.6), radius 0.6*sqrt(2), and order parameter 0
49  GrainSet<dim> test_grain_set_0_new;
50  {
51  test_grain_set_0_new.setGrainIndex(5);
52  test_grain_set_0_new.setOrderParameterIndex(0);
53 
54  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
55  {dealii::Point<dim> p(0.0,0.0); vertex_set[0] = p;}
56  {dealii::Point<dim> p(1.2,0.0); vertex_set[1] = p;}
57  {dealii::Point<dim> p(0.0,1.2); vertex_set[2] = p;}
58  {dealii::Point<dim> p(1.2,1.2); vertex_set[3] = p;}
59 
60  test_grain_set_0_new.addVertexList(vertex_set);
61  }
62 
63  SimplifiedGrainRepresentation<dim> simplified_grain_representation_0_new(test_grain_set_0_new);
64 
65  // Create grain 1_new with center (4.6, 0.6), radius 0.6*sqrt(2), and order parameter 1
66  GrainSet<dim> test_grain_set_1_new;
67  {
68  test_grain_set_1.setGrainIndex(3);
69  test_grain_set_1.setOrderParameterIndex(1);
70 
71  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
72  {dealii::Point<dim> p(4.0,0.0); vertex_set[0] = p;}
73  {dealii::Point<dim> p(5.2,0.0); vertex_set[1] = p;}
74  {dealii::Point<dim> p(4.0,1.2); vertex_set[2] = p;}
75  {dealii::Point<dim> p(5.2,1.2); vertex_set[3] = p;}
76 
77  test_grain_set_1_new.addVertexList(vertex_set);
78  }
79 
80  SimplifiedGrainRepresentation<dim> simplified_grain_representation_1_new(test_grain_set_1_new);
81 
82  // Now build of the vectors of these grain representations
83  std::vector<SimplifiedGrainRepresentation<dim>> old_grain_representations, new_grain_representations;
84  old_grain_representations.push_back(simplified_grain_representation_0);
85  old_grain_representations.push_back(simplified_grain_representation_1);
86  new_grain_representations.push_back(simplified_grain_representation_0_new);
87  new_grain_representations.push_back(simplified_grain_representation_1_new);
88 
89  // Now run the actual test
90  std::cout << "Grain ids before transfer: " << new_grain_representations.at(0).getGrainId() << " " << new_grain_representations.at(1).getGrainId() << std::endl;
91 
92  SimplifiedGrainManipulator<dim> simplified_grain_manipulator;
93  simplified_grain_manipulator.transferGrainIds(old_grain_representations, new_grain_representations);
94 
95  std::cout << "Grain ids after transfer: " << new_grain_representations.at(0).getGrainId() << " " << new_grain_representations.at(1).getGrainId() << std::endl;
96 
97  bool result = false;
98  if ( new_grain_representations.at(0).getGrainId() == 0 and new_grain_representations.at(1).getGrainId() == 1){
99  result = true;
100  }
101  pass = pass & result;
102  }
103 
104  sprintf (buffer, "Test result for 'SimplifiedGrainManipulator::transferGrainIds': %u\n", pass);
105  std::cout << buffer;
106 
107  return pass;
108 }
109 
110 template <int dim,typename T>
112 
113  char buffer[100];
114 
115  std::cout << "\nTesting 'SimplifiedGrainManipulator::reassignGrains'... " << std::endl;
116 
117  bool pass = true;
118 
119  // Subtest 1 (three grains, two order parameters, one grain needs reassignment)
120  {
121  // Create grain 0 with center (0.5, 0.5), radius 0.5*sqrt(2), and order parameter 0
122  GrainSet<dim> test_grain_set_0;
123  {
124  test_grain_set_0.setGrainIndex(0);
125  test_grain_set_0.setOrderParameterIndex(0);
126 
127  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
128  {dealii::Point<dim> p(0.0,0.0); vertex_set[0] = p;}
129  {dealii::Point<dim> p(1.0,0.0); vertex_set[1] = p;}
130  {dealii::Point<dim> p(0.0,1.0); vertex_set[2] = p;}
131  {dealii::Point<dim> p(1.0,1.0); vertex_set[3] = p;}
132 
133  test_grain_set_0.addVertexList(vertex_set);
134  }
135 
136  SimplifiedGrainRepresentation<dim> simplified_grain_representation_0(test_grain_set_0);
137 
138  // Create grain 1 with center (4.5, 0.5), radius 0.5*sqrt(2), and order parameter 1
139  GrainSet<dim> test_grain_set_1;
140  {
141  test_grain_set_1.setGrainIndex(1);
142  test_grain_set_1.setOrderParameterIndex(1);
143 
144  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
145  {dealii::Point<dim> p(4.0,0.0); vertex_set[0] = p;}
146  {dealii::Point<dim> p(5.0,0.0); vertex_set[1] = p;}
147  {dealii::Point<dim> p(4.0,1.0); vertex_set[2] = p;}
148  {dealii::Point<dim> p(5.0,1.0); vertex_set[3] = p;}
149 
150  test_grain_set_1.addVertexList(vertex_set);
151  }
152 
153  SimplifiedGrainRepresentation<dim> simplified_grain_representation_1(test_grain_set_1);
154 
155  // Create grain 2 with center (0.6, 0.6), radius 0.6*sqrt(2), and order parameter 0
156  GrainSet<dim> test_grain_set_2;
157  {
158  test_grain_set_2.setGrainIndex(2);
159  test_grain_set_2.setOrderParameterIndex(0);
160 
161  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
162  {dealii::Point<dim> p(0.0,0.0); vertex_set[0] = p;}
163  {dealii::Point<dim> p(1.2,0.0); vertex_set[1] = p;}
164  {dealii::Point<dim> p(0.0,1.2); vertex_set[2] = p;}
165  {dealii::Point<dim> p(1.2,1.2); vertex_set[3] = p;}
166 
167  test_grain_set_2.addVertexList(vertex_set);
168  }
169 
170  SimplifiedGrainRepresentation<dim> simplified_grain_representation_2(test_grain_set_2);
171 
172 
173  // Now build of the vectors of these grain representations
174  std::vector<SimplifiedGrainRepresentation<dim>> grain_representations;
175  grain_representations.push_back(simplified_grain_representation_0);
176  grain_representations.push_back(simplified_grain_representation_1);
177  grain_representations.push_back(simplified_grain_representation_2);
178 
179  std::vector<unsigned int> order_parameter_id_list;
180  order_parameter_id_list.push_back(0);
181  order_parameter_id_list.push_back(1);
182 
183  // Now run the actual test
184  /*
185  for (unsigned int g=0; g < grain_representations.size(); g++){
186  std::cout << "Grain ops before reassignment: " << grain_representations.at(g).getOrderParameterId() << std::endl;
187  }
188  */
189  SimplifiedGrainManipulator<dim> simplified_grain_manipulator;
190  simplified_grain_manipulator.reassignGrains(grain_representations, 0.5, order_parameter_id_list);
191  /*
192  for (unsigned int g=0; g < grain_representations.size(); g++){
193  std::cout << "Grain ops after reassignment: " << grain_representations.at(g).getOrderParameterId() << std::endl;
194  }
195  */
196  bool result = false;
197  if ( grain_representations.at(0).getOrderParameterId() == 1 and grain_representations.at(1).getOrderParameterId() == 1 and grain_representations.at(2).getOrderParameterId() == 0){
198  result = true;
199  }
200 
201  sprintf (buffer, "Subtest 1 result for 'SimplifiedGrainManipulator::reassignGrains': %u\n", result);
202  std::cout << buffer;
203 
204  pass = pass & result;
205  }
206 
207  // Subtest 2 (three grains, two order parameters, no grains need reassignment)
208  {
209  // Create grain 0 with center (0.5, 0.5), radius 0.5*sqrt(2), and order parameter 0
210  GrainSet<dim> test_grain_set_0;
211  {
212  test_grain_set_0.setGrainIndex(0);
213  test_grain_set_0.setOrderParameterIndex(0);
214 
215  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
216  {dealii::Point<dim> p(0.0,0.0); vertex_set[0] = p;}
217  {dealii::Point<dim> p(1.0,0.0); vertex_set[1] = p;}
218  {dealii::Point<dim> p(0.0,1.0); vertex_set[2] = p;}
219  {dealii::Point<dim> p(1.0,1.0); vertex_set[3] = p;}
220 
221  test_grain_set_0.addVertexList(vertex_set);
222  }
223 
224  SimplifiedGrainRepresentation<dim> simplified_grain_representation_0(test_grain_set_0);
225 
226  // Create grain 1 with center (4.5, 0.5), radius 0.5*sqrt(2), and order parameter 1
227  GrainSet<dim> test_grain_set_1;
228  {
229  test_grain_set_1.setGrainIndex(1);
230  test_grain_set_1.setOrderParameterIndex(1);
231 
232  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
233  {dealii::Point<dim> p(4.0,0.0); vertex_set[0] = p;}
234  {dealii::Point<dim> p(5.0,0.0); vertex_set[1] = p;}
235  {dealii::Point<dim> p(4.0,1.0); vertex_set[2] = p;}
236  {dealii::Point<dim> p(5.0,1.0); vertex_set[3] = p;}
237 
238  test_grain_set_1.addVertexList(vertex_set);
239  }
240 
241  SimplifiedGrainRepresentation<dim> simplified_grain_representation_1(test_grain_set_1);
242 
243  // Create grain 2 with center (2.5, 2.5), radius 0.5*sqrt(2), and order parameter 0
244  GrainSet<dim> test_grain_set_2;
245  {
246  test_grain_set_2.setGrainIndex(2);
247  test_grain_set_2.setOrderParameterIndex(0);
248 
249  std::vector<dealii::Point<dim>> vertex_set(dealii::Utilities::fixed_power<dim>(2.0));
250  {dealii::Point<dim> p(2.0,2.0); vertex_set[0] = p;}
251  {dealii::Point<dim> p(3.0,2.0); vertex_set[1] = p;}
252  {dealii::Point<dim> p(2.0,3.0); vertex_set[2] = p;}
253  {dealii::Point<dim> p(3.0,3.0); vertex_set[3] = p;}
254 
255  test_grain_set_2.addVertexList(vertex_set);
256  }
257 
258  SimplifiedGrainRepresentation<dim> simplified_grain_representation_2(test_grain_set_2);
259 
260 
261  // Now build of the vectors of these grain representations
262  std::vector<SimplifiedGrainRepresentation<dim>> grain_representations;
263  grain_representations.push_back(simplified_grain_representation_0);
264  grain_representations.push_back(simplified_grain_representation_1);
265  grain_representations.push_back(simplified_grain_representation_2);
266 
267  std::vector<unsigned int> order_parameter_id_list;
268  order_parameter_id_list.push_back(0);
269  order_parameter_id_list.push_back(1);
270 
271  // Now run the actual test
272  /*
273  for (unsigned int g=0; g < grain_representations.size(); g++){
274  std::cout << "Grain ops before reassignment: " << grain_representations.at(g).getOrderParameterId() << std::endl;
275  }
276  */
277 
278  SimplifiedGrainManipulator<dim> simplified_grain_manipulator;
279  simplified_grain_manipulator.reassignGrains(grain_representations, 0.5, order_parameter_id_list);
280 
281  /*
282  for (unsigned int g=0; g < grain_representations.size(); g++){
283  std::cout << "Grain ops after reassignment: " << grain_representations.at(g).getOrderParameterId() << std::endl;
284  }
285  */
286 
287  bool result = false;
288  if ( grain_representations.at(0).getOrderParameterId() == 0 and grain_representations.at(1).getOrderParameterId() == 1 and grain_representations.at(2).getOrderParameterId() == 0){
289  result = true;
290  }
291 
292  sprintf (buffer, "Subtest 2 result for 'SimplifiedGrainManipulator::reassignGrains': %u\n", result);
293  std::cout << buffer;
294 
295  pass = pass & result;
296  }
297 
298  sprintf (buffer, "Test result for 'SimplifiedGrainManipulator::reassignGrains': %u\n", pass);
299  std::cout << buffer;
300 
301  return pass;
302 }
void reassignGrains(std::vector< SimplifiedGrainRepresentation< dim >> &grain_representations, double buffer_distance, std::vector< unsigned int > order_parameter_id_list)
void setOrderParameterIndex(unsigned int _order_parameter_index)
Definition: FloodFiller.h:38
void transferGrainIds(const std::vector< SimplifiedGrainRepresentation< dim >> &old_grain_representations, std::vector< SimplifiedGrainRepresentation< dim >> &new_grain_representations) const
void addVertexList(std::vector< dealii::Point< dim >> _vertices)
Definition: FloodFiller.h:48
bool test_SimplifiedGrainManipulator_transferGrainIds()
bool test_SimplifiedGrainManipulator_reassignGrains()
void setGrainIndex(unsigned int _grain_index)
Definition: FloodFiller.h:28