MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_orthotropic_material.cpp
Go to the documentation of this file.
1 // Catch2 includes
2 #include "catch.hpp"
3 
4 // libMesh includes
5 #include "libmesh/point.h"
6 
7 // MAST includes
8 #include "base/parameter.h"
11 
12 // Custom includes
13 #include "test_helpers.h"
14 
15 // TODO: Need to test temperature dependent material property
16 // TODO: Need to test stress dependent material property (plasticity)
17 // TODO: Need to test material property where field functions are not constant
18 // TODO: Need to implement tests for 3D elements
19 
20 
21 TEST_CASE("constant_orthotropic_thermoelastic_material_1d",
22  "[orthotropic],[material],[constant],[1D],[thermoelastic]")
23 {
24  MAST::Parameter alpha11("alpha11_param", 5.43e-05); // Coefficient of thermal expansion
25  MAST::ConstantFieldFunction alpha11_f("alpha11_expansion", alpha11);
26 
27  // Initialize the material
29 
30  material.add(alpha11_f);
31 
32  REQUIRE( material.contains("alpha11_expansion") );
33 
34  SECTION("1D material thermal expansion matrix")
35  {
36  const MAST::FieldFunction<RealMatrixX>& mat_texp =
37  material.thermal_expansion_matrix(1);
38 
39  const libMesh::Point point(2.3, 3.1, 5.2);
40  const Real time = 2.34;
41  RealMatrixX D_texp;
42  mat_texp(point, time, D_texp);
43 
44  // Hard-coded in the true value of material thermal expansion matrix
45  RealMatrixX D_texp_true = RealMatrixX::Zero(2,1);
46  D_texp_true(0,0) = 5.43e-05;
47 
48  // Convert the test and truth Eigen::Matrix objects to std::vector
49  // since Catch2 has built in methods to compare vectors
50  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_texp);
51  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_texp_true);
52 
53  // Floating point approximations are diffcult to compare since the
54  // values typically aren't exactly equal due to numerical error.
55  // Therefore, we use the Approx comparison instead of Equals
56  CHECK_THAT( test, Catch::Approx<double>(truth) );
57  }
58 }
59 
60 
61 TEST_CASE("constant_orthotropic_structural_material_1d",
62  "[orthotropic],[material],[constant],[1D],[structural]")
63 {
64  // Define Material Properties as MAST Parameters
65  MAST::Parameter E11("E11_param", 72.0e9); // Modulus of Elasticity
66  MAST::Parameter G12("G12_param", 68.3e9); // Shear Modulus
67 
68  // Create field functions to dsitribute these constant parameters throughout the model
69  MAST::ConstantFieldFunction E11_f("E11", E11);
70  MAST::ConstantFieldFunction G12_f("G12", G12);
71 
72  // Initialize the material
74 
75  // Add the material property constant field functions to the material card
76  material.add(E11_f);
77  material.add(G12_f);
78 
79 
80  // Ensure that the material properties were added before doing other checks
81  REQUIRE( material.contains("E11") );
82  REQUIRE( material.contains("G12") );
83 
84  SECTION("1D material constitutive matrix")
85  {
86  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
87  material.stiffness_matrix(1);
88 
89  const libMesh::Point point(2.3, 3.1, 5.2);
90  const Real time = 2.34;
91  RealMatrixX D;
92  mat_stiffness(point, time, D);
93 
94  // Hard-code in the true value of the material stiffness
95  RealMatrixX D_true = RealMatrixX::Zero(2,2);
96  D_true(0,0) = 72.0e9;
97  D_true(1,1) = 68.3e9;
98 
99  // Convert the test and truth Eigen::Matrix objects to std::vector
100  // since Catch2 has built in methods to compare vectors
101  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D);
102  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_true);
103 
104  // Floating point approximations are diffcult to compare since the
105  // values typically aren't exactly equal due to numerical error.
106  // Therefore, we use the Approx comparison instead of Equals
107  CHECK_THAT( test, Catch::Approx<double>(truth) );
108  }
109 }
110 
111 
112 TEST_CASE("constant_orthotropic_heat_transfer_material_1d",
113  "[orthotropic],[material],[1D],[heat_transfer][constant]")
114 {
115  // Define Material Properties as MAST Parameters
116  MAST::Parameter k11("k11_param", 237.0); // Thermal Conductivity
117 
118  // Create field functions to dsitribute these constant parameters throughout the model
119  MAST::ConstantFieldFunction k11_f("k11_th", k11);
120 
121  // Initialize the material
123 
124  // Add the material property constant field functions to the material card
125  material.add(k11_f);
126 
127  // Ensure that the material properties were added before doing other checks
128  REQUIRE( material.contains("k11_th") );
129 
130  SECTION("material depends on the parameters that it should")
131  {
132  CHECK( material.depends_on(k11) );
133  }
134 
135  SECTION("material does not depend on other parameters")
136  {
137  MAST::Parameter dummy("dummy", 1.0);
138  CHECK_FALSE( material.depends_on(dummy) );
139  }
140 
141  SECTION("1D thermal conductivity matrix")
142  {
143  const MAST::FieldFunction<RealMatrixX>& mat_conduct =
144  material.conductance_matrix(1);
145 
146  const libMesh::Point point(2.3, 3.1, 5.2);
147  const Real time = 2.34;
148  RealMatrixX D_k;
149  mat_conduct(point, time, D_k);
150 
151  // Hard-coded values for thermal conductivity matrix
152  RealMatrixX D_k_true = RealMatrixX::Identity(1,1);
153  D_k_true *= 237.0;
154 
155  // Convert the test and truth Eigen::Matrix objects to std::vector
156  // since Catch2 has built in methods to compare vectors
157  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_k);
158  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_k_true);
159 
160  // Floating point approximations are diffcult to compare since the
161  // values typically aren't exactly equal due to numerical error.
162  // Therefore, we use the Approx comparison instead of Equals
163  CHECK_THAT( test, Catch::Approx<double>(truth) );
164  }
165 }
166 
167 
168 TEST_CASE("constant_orthotropic_transient_heat_transfer_material_1d",
169  "[orthotropic],[material],[1D],[heat_transfer],[constant],[transient]")
170 {
171  // Define Material Properties as MAST Parameters
172  MAST::Parameter rho("rho_param", 1234.5); // Density
173  MAST::Parameter cp("cp_param", 908.0); // Specific Heat Capacity
174 
175  // Create field functions to dsitribute these constant parameters throughout the model
176  MAST::ConstantFieldFunction rho_f("rho", rho);
177  MAST::ConstantFieldFunction cp_f("cp", cp);
178 
179  // Initialize the material
181 
182  // Add the material property constant field functions to the material card
183  material.add(rho_f);
184  material.add(cp_f);
185 
186  // Ensure that the material properties were added before doing other checks
187  REQUIRE( material.contains("rho") );
188  REQUIRE( material.contains("cp") );
189 
190  SECTION("material depends on the parameters that it should")
191  {
192  CHECK( material.depends_on(rho) );
193  CHECK( material.depends_on(cp) );
194  }
195 
196  SECTION("material does not depend on other parameters")
197  {
198  MAST::Parameter dummy("dummy", 1.0);
199  CHECK_FALSE( material.depends_on(dummy) );
200  }
201 
202  SECTION("1D thermal capacitance matrix")
203  {
204  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
205  material.capacitance_matrix(1);
206 
207  const libMesh::Point point(2.3, 3.1, 5.2);
208  const Real time = 2.34;
209  RealMatrixX D_cp;
210  mat_capacit(point, time, D_cp);
211 
212  // Hard-coded values for thermal conductivity matrix
213  RealMatrixX D_cp_true = RealMatrixX::Identity(1,1);
214  D_cp_true *= (908.0 * 1234.5);
215 
216  // Convert the test and truth Eigen::Matrix objects to std::vector
217  // since Catch2 has built in methods to compare vectors
218  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_cp);
219  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_cp_true);
220 
221  // Floating point approximations are diffcult to compare since the
222  // values typically aren't exactly equal due to numerical error.
223  // Therefore, we use the Approx comparison instead of Equals
224  CHECK_THAT( test, Catch::Approx<double>(truth) );
225  }
226 }
227 
228 
229 TEST_CASE("constant_orthotropic_thermoelastic_material_2d",
230  "[orthotropic],[material],[2D],[thermoelastic][constant]")
231 {
232  MAST::Parameter alpha11("alpha11_param", 5.43e-05); // Coefficient of thermal expansion
233  MAST::Parameter alpha22("alpha22_param", 8.79e-06); // Coefficient of thermal expansion
234 
235  MAST::ConstantFieldFunction alpha11_f("alpha11_expansion", alpha11);
236  MAST::ConstantFieldFunction alpha22_f("alpha22_expansion", alpha22);
237 
238  // Initialize the material
240 
241  material.add(alpha11_f);
242  material.add(alpha22_f);
243 
244  REQUIRE( material.contains("alpha11_expansion") );
245  REQUIRE( material.contains("alpha22_expansion") );
246 
247  SECTION("material does not depend on other parameters")
248  {
249  MAST::Parameter dummy("dummy", 1.0);
250  CHECK_FALSE( material.depends_on(dummy) );
251  }
252 
253  SECTION("2D material thermal expansion matrix")
254  {
255  const MAST::FieldFunction<RealMatrixX>& mat_texp =
256  material.thermal_expansion_matrix(2);
257 
258  const libMesh::Point point(2.3, 3.1, 5.2);
259  const Real time = 2.34;
260  RealMatrixX D_texp;
261  mat_texp(point, time, D_texp);
262 
263  // Hard-coded in the true value of material thermal expansion matrix
264  RealMatrixX D_texp_true = RealMatrixX::Zero(3,1);
265  D_texp_true(0,0) = 5.43e-05;
266  D_texp_true(1,0) = 8.79e-06;
267 
268  // Convert the test and truth Eigen::Matrix objects to std::vector
269  // since Catch2 has built in methods to compare vectors
270  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_texp);
271  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_texp_true);
272 
273  // Floating point approximations are diffcult to compare since the
274  // values typically aren't exactly equal due to numerical error.
275  // Therefore, we use the Approx comparison instead of Equals
276  CHECK_THAT( test, Catch::Approx<double>(truth) );
277  }
278 }
279 
280 
281 TEST_CASE("constant_orthotropic_structural_material_2d",
282  "[orthotropic],[material],[2D],[structural],[constant]")
283 {
284  // Define Material Properties as MAST Parameters
285  MAST::Parameter E11("E11", 72.0e9); // Modulus of Elasticity
286  MAST::Parameter E22("E22", 55.0e9);
287  MAST::Parameter E33("E33", 33.4e9);
288  MAST::Parameter nu12("nu12", 0.33); // Poisson's ratio
289  MAST::Parameter nu23("nu23", 0.25);
290  MAST::Parameter nu31("nu31", 0.45);
291  MAST::Parameter G12("G12", 68.3e9);
292 
293  // Create field functions to dsitribute these constant parameters throughout the model
294  MAST::ConstantFieldFunction E11_f("E11", E11);
295  MAST::ConstantFieldFunction E22_f("E22", E22);
296  MAST::ConstantFieldFunction E33_f("E33", E33);
297  MAST::ConstantFieldFunction nu12_f("nu12", nu12);
298  MAST::ConstantFieldFunction nu23_f("nu23", nu23);
299  MAST::ConstantFieldFunction nu31_f("nu31", nu31);
300  MAST::ConstantFieldFunction G12_f("G12", G12);
301 
302  // Initialize the material
304 
305  // Add the material property constant field functions to the material card
306  material.add(E11_f);
307  material.add(E22_f);
308  material.add(E33_f);
309  material.add(nu12_f);
310  material.add(nu23_f);
311  material.add(nu31_f);
312  material.add(G12_f);
313 
314  // Ensure that the material properties were added before doing other checks
315  REQUIRE( material.contains("E11") );
316  REQUIRE( material.contains("E22") );
317  REQUIRE( material.contains("E33") );
318  REQUIRE( material.contains("nu12") );
319  REQUIRE( material.contains("nu23") );
320  REQUIRE( material.contains("nu31") );
321  REQUIRE( material.contains("G12") );
322 
323  SECTION("material depends on the parameters that it should")
324  {
325  CHECK( material.depends_on(E11) );
326  CHECK( material.depends_on(E22) );
327  CHECK( material.depends_on(E33) );
328  CHECK( material.depends_on(nu12) );
329  CHECK( material.depends_on(nu23) );
330  CHECK( material.depends_on(nu31) );
331  CHECK( material.depends_on(G12) );
332  }
333 
334  SECTION("material does not depend on other parameters")
335  {
336  MAST::Parameter dummy("dummy", 1.0);
337  CHECK_FALSE( material.depends_on(dummy) );
338  }
339 
340  SECTION("2D plane stress material constitutive matrix")
341  {
342  const bool is_plane_stress = true;
343  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
344  material.stiffness_matrix(2, is_plane_stress);
345 
346  const libMesh::Point point(2.3, 3.1, 5.2);
347  const Real time = 2.34;
348  RealMatrixX K_mat;
349  mat_stiffness(point, time, K_mat);
350 
351  // Check for symmetry
352 
353  // Convert the test and truth Eigen::Matrix objects to std::vector
354  // since Catch2 has built in methods to compare vectors
355  std::vector<double> test = TEST::eigen_matrix_to_std_vector(K_mat.transpose());
356  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(K_mat);
357 
358  // Floating point approximations are diffcult to compare since the
359  // values typically aren't exactly equal due to numerical error.
360  // Therefore, we use the Approx comparison instead of Equals
361  REQUIRE_THAT( test, Catch::Approx<double>(truth) );
362 
363  // Check for positive definiteness
364  SelfAdjointEigenSolver<RealMatrixX> eigensolver(K_mat);
365  if (eigensolver.info() != Success)
366  {
367  std::cout << "eigensolver failed to converge!" << std::endl;
368  REQUIRE(false); // Force test failure.
369  }
370  RealVectorX eigenvalues = eigensolver.eigenvalues();
371  for (uint i=0; i<K_mat.cols(); i++)
372  {
373  REQUIRE( eigenvalues(i) > 0.0 );
374  }
375 
376  // Hard-code in the true value of the material stiffness
377  RealMatrixX K_mat_true = RealMatrixX::Zero(3,3);
378  K_mat_true(0,0) = 7.853296066534869e+10;
379  K_mat_true(1,1) = 5.999045606380803e+10;
380  K_mat_true(2,2) = 6.830000000000000e+10;
381  K_mat_true(0,1) = K_mat_true(1,0) =1.979685050105665e+10;
382 
383  // Convert the test and truth Eigen::Matrix objects to std::vector
384  // since Catch2 has built in methods to compare vectors
385  test = TEST::eigen_matrix_to_std_vector(K_mat);
386  truth = TEST::eigen_matrix_to_std_vector(K_mat_true);
387 
388  // Floating point approximations are diffcult to compare since the
389  // values typically aren't exactly equal due to numerical error.
390  // Therefore, we use the Approx comparison instead of Equals
391  CHECK_THAT( test, Catch::Approx<double>(truth) );
392  }
393 
394  SECTION("2D plane strain material constitutive matrix")
395  {
396  const bool is_plane_stress = false;
397  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
398  material.stiffness_matrix(2, is_plane_stress);
399 
400  const libMesh::Point point(2.3, 3.1, 5.2);
401  const Real time = 2.34;
402  RealMatrixX K_mat;
403  mat_stiffness(point, time, K_mat);
404 
405  // Check for symmetry
406 
407  // Convert the test and truth Eigen::Matrix objects to std::vector
408  // since Catch2 has built in methods to compare vectors
409  std::vector<double> test = TEST::eigen_matrix_to_std_vector(K_mat.transpose());
410  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(K_mat);
411 
412  // Floating point approximations are diffcult to compare since the
413  // values typically aren't exactly equal due to numerical error.
414  // Therefore, we use the Approx comparison instead of Equals
415  REQUIRE_THAT( test, Catch::Approx<double>(truth) );
416 
417  // Check for positive definiteness
418  SelfAdjointEigenSolver<RealMatrixX> eigensolver(K_mat);
419  if (eigensolver.info() != Success)
420  {
421  std::cout << "eigensolver failed to converge!" << std::endl;
422  REQUIRE(false); // Force test failure.
423  }
424  RealVectorX eigenvalues = eigensolver.eigenvalues();
425  for (uint i=0; i<K_mat.cols(); i++)
426  {
427  REQUIRE( eigenvalues(i) > 0.0 );
428  }
429 
430  // Hard-code in the true value of the material stiffness
431  RealMatrixX K_mat_true = RealMatrixX::Zero(3,3);
432  K_mat_true(0,0) = 1.881848591463047e11;
433  K_mat_true(1,1) = 0.841961884847417e11;
434  K_mat_true(2,2) = 0.683000000000000e11;
435  K_mat_true(0,1) = K_mat_true(1,0) = 0.713158228712175e11;
436 
437  // Convert the test and truth Eigen::Matrix objects to std::vector
438  // since Catch2 has built in methods to compare vectors
439  test = TEST::eigen_matrix_to_std_vector(K_mat);
440  truth = TEST::eigen_matrix_to_std_vector(K_mat_true);
441 
442  // Floating point approximations are diffcult to compare since the
443  // values typically aren't exactly equal due to numerical error.
444  // Therefore, we use the Approx comparison instead of Equals
445  CHECK_THAT( test, Catch::Approx<double>(truth) );
446  }
447 
448  SECTION("transverse shear material stiffness matrix")
449  {
450  const MAST::FieldFunction<RealMatrixX>& mat_trans_shear =
452 
453  const libMesh::Point point(2.3, 3.1, 5.2);
454  const Real time = 2.34;
455  RealMatrixX D_trans_shear;
456  mat_trans_shear(point, time, D_trans_shear);
457 
458  // Hard-coded in the true value of the material transverse shear stiffness
459  RealMatrixX D_trans_shear_true = RealMatrixX::Zero(2,2);
460  D_trans_shear_true(0,0) = 6.830000000000000e+10;
461  D_trans_shear_true(1,1) = 6.830000000000000e+10;
462 
463  // Convert the test and truth Eigen::Matrix objects to std::vector
464  // since Catch2 has built in methods to compare vectors
465  std::vector<double> test(D_trans_shear.data(), D_trans_shear.data()+D_trans_shear.rows()*D_trans_shear.cols());
466  std::vector<double> truth(D_trans_shear_true.data(), D_trans_shear_true.data()+D_trans_shear_true.rows()*D_trans_shear_true.cols());
467 
468  // Floating point approximations are diffcult to compare since the
469  // values typically aren't exactly equal due to numerical error.
470  // Therefore, we use the Approx comparison instead of Equals
471  CHECK_THAT( test, Catch::Approx<double>(truth) );
472  }
473 }
474 
475 
476 TEST_CASE("constant_orthotropic_heat_transfer_material_2d",
477  "[orthotropic],[material],[2D],[heat_transfer][constant]")
478 {
479  // Define Material Properties as MAST Parameters
480  MAST::Parameter k11("k11_param", 237.0); // Thermal Conductivity
481  MAST::Parameter k22("k22_param", 542.0); // Thermal Conductivity
482 
483  // Create field functions to dsitribute these constant parameters throughout the model
484  MAST::ConstantFieldFunction k11_f("k11_th", k11);
485  MAST::ConstantFieldFunction k22_f("k22_th", k22);
486 
487  // Initialize the material
489 
490  // Add the material property constant field functions to the material card
491  material.add(k11_f);
492  material.add(k22_f);
493 
494  // Ensure that the material properties were added before doing other checks
495  REQUIRE( material.contains("k11_th") );
496  REQUIRE( material.contains("k22_th") );
497 
498  SECTION("material depends on the parameters that it should")
499  {
500  CHECK( material.depends_on(k11) );
501  CHECK( material.depends_on(k22) );
502  }
503 
504  SECTION("material does not depend on other parameters")
505  {
506  MAST::Parameter dummy("dummy", 1.0);
507  CHECK_FALSE( material.depends_on(dummy) );
508  }
509 
510  SECTION("2D thermal conductivity matrix")
511  {
512  const MAST::FieldFunction<RealMatrixX>& mat_conduct =
513  material.conductance_matrix(2);
514 
515  const libMesh::Point point(2.3, 3.1, 5.2);
516  const Real time = 2.34;
517  RealMatrixX D_k;
518  mat_conduct(point, time, D_k);
519 
520  // Hard-coded values for thermal conductivity matrix
521  RealMatrixX D_k_true = RealMatrixX::Identity(2,2);
522  D_k_true(0,0) = 237.0;
523  D_k_true(1,1) = 542.0;
524 
525  // Convert the test and truth Eigen::Matrix objects to std::vector
526  // since Catch2 has built in methods to compare vectors
527  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_k);
528  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_k_true);
529 
530  // Floating point approximations are diffcult to compare since the
531  // values typically aren't exactly equal due to numerical error.
532  // Therefore, we use the Approx comparison instead of Equals
533  CHECK_THAT( test, Catch::Approx<double>(truth) );
534  }
535 }
536 
537 
538 TEST_CASE("constant_orthotropic_transient_heat_transfer_material_2d",
539  "[orthotropic],[material],[2D],[heat_transfer],[constant],[transient]")
540 {
541  // Define Material Properties as MAST Parameters
542  MAST::Parameter rho("rho_param", 1234.5); // Density
543  MAST::Parameter cp("cp_param", 908.0); // Specific Heat Capacity
544 
545  // Create field functions to dsitribute these constant parameters throughout the model
546  MAST::ConstantFieldFunction rho_f("rho", rho);
547  MAST::ConstantFieldFunction cp_f("cp", cp);
548 
549  // Initialize the material
551 
552  // Add the material property constant field functions to the material card
553  material.add(rho_f);
554  material.add(cp_f);
555 
556  // Ensure that the material properties were added before doing other checks
557  REQUIRE( material.contains("rho") );
558  REQUIRE( material.contains("cp") );
559 
560  SECTION("material depends on the parameters that it should")
561  {
562  CHECK( material.depends_on(rho) );
563  CHECK( material.depends_on(cp) );
564  }
565 
566  SECTION("material does not depend on other parameters")
567  {
568  MAST::Parameter dummy("dummy", 1.0);
569  CHECK_FALSE( material.depends_on(dummy) );
570  }
571 
572  SECTION("2D thermal capacitance matrix")
573  {
574  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
575  material.capacitance_matrix(2);
576 
577  const libMesh::Point point(2.3, 3.1, 5.2);
578  const Real time = 2.34;
579  RealMatrixX D_cp;
580  mat_capacit(point, time, D_cp);
581 
582  // Hard-coded values for thermal conductivity matrix
583  RealMatrixX D_cp_true = RealMatrixX::Identity(1,1);
584  D_cp_true *= (908.0 * 1234.5);
585 
586  // Convert the test and truth Eigen::Matrix objects to std::vector
587  // since Catch2 has built in methods to compare vectors
588  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_cp);
589  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_cp_true);
590 
591  // Floating point approximations are diffcult to compare since the
592  // values typically aren't exactly equal due to numerical error.
593  // Therefore, we use the Approx comparison instead of Equals
594  CHECK_THAT( test, Catch::Approx<double>(truth) );
595  }
596 }
597 
598 
599 TEST_CASE("constant_orthotropic_thermoelastic_material_3d",
600  "[orthotropic],[material],[3D],[thermoelastic][constant]")
601 {
602  MAST::Parameter alpha11("alpha11_param", 5.43e-05);
603  MAST::Parameter alpha22("alpha22_param", 8.79e-06);
604  MAST::Parameter alpha33("alpha33_param", 2.14e-05);
605 
606  MAST::ConstantFieldFunction alpha11_f("alpha11_expansion", alpha11);
607  MAST::ConstantFieldFunction alpha22_f("alpha22_expansion", alpha22);
608  MAST::ConstantFieldFunction alpha33_f("alpha33_expansion", alpha33);
609 
610  // Initialize the material
612 
613  material.add(alpha11_f);
614  material.add(alpha22_f);
615  material.add(alpha33_f);
616 
617  REQUIRE( material.contains("alpha11_expansion") );
618  REQUIRE( material.contains("alpha22_expansion") );
619  REQUIRE( material.contains("alpha33_expansion") );
620 
621  SECTION("material does not depend on other parameters")
622  {
623  MAST::Parameter dummy("dummy", 1.0);
624  CHECK_FALSE( material.depends_on(dummy) );
625  }
626 
627  SECTION("3D material thermal expansion matrix")
628  {
629  const MAST::FieldFunction<RealMatrixX>& mat_texp =
630  material.thermal_expansion_matrix(3);
631 
632  const libMesh::Point point(2.3, 3.1, 5.2);
633  const Real time = 2.34;
634  RealMatrixX D_texp;
635  mat_texp(point, time, D_texp);
636 
637  // Hard-coded in the true value of material thermal expansion matrix
638  RealMatrixX D_texp_true = RealMatrixX::Zero(6,1);
639  D_texp_true(0,0) = 5.43e-05;
640  D_texp_true(1,0) = 8.79e-06;
641  D_texp_true(2,0) = 2.14e-05;
642 
643  // Convert the test and truth Eigen::Matrix objects to std::vector
644  // since Catch2 has built in methods to compare vectors
645  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_texp);
646  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_texp_true);
647 
648  // Floating point approximations are diffcult to compare since the
649  // values typically aren't exactly equal due to numerical error.
650  // Therefore, we use the Approx comparison instead of Equals
651  CHECK_THAT( test, Catch::Approx<double>(truth) );
652  }
653 }
654 
655 
656 TEST_CASE("constant_orthotropic_structural_material_3d",
657  "[orthotropic],[material],[3D],[structural],[constant]")
658 {
659  // Define Material Properties as MAST Parameters
660  MAST::Parameter E11("E11", 72.0e9); // Modulus of Elasticity
661  MAST::Parameter E22("E22", 55.0e9);
662  MAST::Parameter E33("E33", 33.4e9);
663  MAST::Parameter nu12("nu12", 0.33); // Poisson's ratio
664  MAST::Parameter nu23("nu23", 0.25);
665  MAST::Parameter nu31("nu31", 0.45);
666  MAST::Parameter G12("G12", 68.3e9); // Shear Modulus
667  MAST::Parameter G23("G23", 74.5e9);
668  MAST::Parameter G13("G13", 55.1e9);
669 
670  // Create field functions to dsitribute these constant parameters throughout the model
671  MAST::ConstantFieldFunction E11_f("E11", E11);
672  MAST::ConstantFieldFunction E22_f("E22", E22);
673  MAST::ConstantFieldFunction E33_f("E33", E33);
674  MAST::ConstantFieldFunction nu12_f("nu12", nu12);
675  MAST::ConstantFieldFunction nu23_f("nu23", nu23);
676  MAST::ConstantFieldFunction nu31_f("nu31", nu31);
677  MAST::ConstantFieldFunction G12_f("G12", G12);
678  MAST::ConstantFieldFunction G23_f("G23", G23);
679  MAST::ConstantFieldFunction G13_f("G13", G13);
680 
681  // Initialize the material
683 
684  // Add the material property constant field functions to the material card
685  material.add(E11_f);
686  material.add(E22_f);
687  material.add(E33_f);
688  material.add(nu12_f);
689  material.add(nu23_f);
690  material.add(nu31_f);
691  material.add(G12_f);
692  material.add(G23_f);
693  material.add(G13_f);
694 
695  // Ensure that the material properties were added before doing other checks
696  REQUIRE( material.contains("E11") );
697  REQUIRE( material.contains("E22") );
698  REQUIRE( material.contains("E33") );
699  REQUIRE( material.contains("nu12") );
700  REQUIRE( material.contains("nu23") );
701  REQUIRE( material.contains("nu31") );
702  REQUIRE( material.contains("G12") );
703  REQUIRE( material.contains("G23") );
704  REQUIRE( material.contains("G13") );
705 
706  SECTION("material depends on the parameters that it should")
707  {
708  CHECK( material.depends_on(E11) );
709  CHECK( material.depends_on(E22) );
710  CHECK( material.depends_on(E33) );
711  CHECK( material.depends_on(nu12) );
712  CHECK( material.depends_on(nu23) );
713  CHECK( material.depends_on(nu31) );
714  CHECK( material.depends_on(G12) );
715  CHECK( material.depends_on(G23) );
716  CHECK( material.depends_on(G13) );
717  }
718 
719  SECTION("material does not depend on other parameters")
720  {
721  MAST::Parameter dummy("dummy", 1.0);
722  CHECK_FALSE( material.depends_on(dummy) );
723  }
724 
725  SECTION("3D material constitutive matrix")
726  {
727  const bool is_plane_stress = true;
728  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
729  material.stiffness_matrix(3);
730 
731  const libMesh::Point point(2.3, 3.1, 5.2);
732  const Real time = 2.34;
733  RealMatrixX K_mat;
734  mat_stiffness(point, time, K_mat);
735 
736  // Check for symmetry
737 
738  // Convert the test and truth Eigen::Matrix objects to std::vector
739  // since Catch2 has built in methods to compare vectors
740  std::vector<double> test = TEST::eigen_matrix_to_std_vector(K_mat.transpose());
741  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(K_mat);
742 
743  // Floating point approximations are diffcult to compare since the
744  // values typically aren't exactly equal due to numerical error.
745  // Therefore, we use the Approx comparison instead of Equals
746  REQUIRE_THAT( test, Catch::Approx<double>(truth) );
747 
748  // Check for positive definiteness
749  SelfAdjointEigenSolver<RealMatrixX> eigensolver(K_mat);
750  if (eigensolver.info() != Success)
751  {
752  std::cout << "eigensolver failed to converge!" << std::endl;
753  REQUIRE(false); // Force test failure.
754  }
755  RealVectorX eigenvalues = eigensolver.eigenvalues();
756  for (uint i=0; i<K_mat.cols(); i++)
757  {
758  REQUIRE( eigenvalues(i) > 0.0 );
759  }
760 
761  // Hard-code in the true value of the material stiffness
762  RealMatrixX K_mat_true = RealMatrixX::Zero(6,6);
763  K_mat_true(0,0) = 1.881848591463047e+11;
764  K_mat_true(1,1) = 0.841961884847417e+11;
765  K_mat_true(2,2) = 0.831923864531179e+11;
766  K_mat_true(3,3) = 0.683000000000000e+11;
767  K_mat_true(4,4) = 0.745000000000000e+11;
768  K_mat_true(5,5) = 0.551000000000000e+11;
769 
770  K_mat_true(0,1) = K_mat_true(1,0) = 0.713158228712175e+11;
771  K_mat_true(0,2) = K_mat_true(2,0) = 0.955102251790129e+11;
772  K_mat_true(1,2) = K_mat_true(2,1) = 0.448746325438223e+11;
773 
774  // Convert the test and truth Eigen::Matrix objects to std::vector
775  // since Catch2 has built in methods to compare vectors
776  test = TEST::eigen_matrix_to_std_vector(K_mat);
777  truth = TEST::eigen_matrix_to_std_vector(K_mat_true);
778 
779  // Floating point approximations are diffcult to compare since the
780  // values typically aren't exactly equal due to numerical error.
781  // Therefore, we use the Approx comparison instead of Equals
782  CHECK_THAT( test, Catch::Approx<double>(truth) );
783  }
784 }
785 
786 
787 TEST_CASE("constant_orthotropic_heat_transfer_material_3d",
788  "[orthotropic],[material],[3D],[heat_transfer][constant]")
789 {
790  // Define Material Properties as MAST Parameters
791  MAST::Parameter k11("k11_param", 237.0); // Thermal Conductivity
792  MAST::Parameter k22("k22_param", 542.0); // Thermal Conductivity
793  MAST::Parameter k33("k33_param", 444.4); // Thermal Conductivity
794 
795  // Create field functions to dsitribute these constant parameters throughout the model
796  MAST::ConstantFieldFunction k11_f("k11_th", k11);
797  MAST::ConstantFieldFunction k22_f("k22_th", k22);
798  MAST::ConstantFieldFunction k33_f("k33_th", k33);
799 
800  // Initialize the material
802 
803  // Add the material property constant field functions to the material card
804  material.add(k11_f);
805  material.add(k22_f);
806  material.add(k33_f);
807 
808  // Ensure that the material properties were added before doing other checks
809  REQUIRE( material.contains("k11_th") );
810  REQUIRE( material.contains("k22_th") );
811  REQUIRE( material.contains("k33_th") );
812 
813  SECTION("material depends on the parameters that it should")
814  {
815  CHECK( material.depends_on(k11) );
816  CHECK( material.depends_on(k22) );
817  CHECK( material.depends_on(k33) );
818  }
819 
820  SECTION("material does not depend on other parameters")
821  {
822  MAST::Parameter dummy("dummy", 1.0);
823  CHECK_FALSE( material.depends_on(dummy) );
824  }
825 
826  SECTION("3D thermal conductivity matrix")
827  {
828  const MAST::FieldFunction<RealMatrixX>& mat_conduct =
829  material.conductance_matrix(3);
830 
831  const libMesh::Point point(2.3, 3.1, 5.2);
832  const Real time = 2.34;
833  RealMatrixX D_k;
834  mat_conduct(point, time, D_k);
835 
836  // Hard-coded values for thermal conductivity matrix
837  RealMatrixX D_k_true = RealMatrixX::Identity(3,3);
838  D_k_true(0,0) = 237.0;
839  D_k_true(1,1) = 542.0;
840  D_k_true(2,2) = 444.4;
841 
842  // Convert the test and truth Eigen::Matrix objects to std::vector
843  // since Catch2 has built in methods to compare vectors
844  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_k);
845  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_k_true);
846 
847  // Floating point approximations are diffcult to compare since the
848  // values typically aren't exactly equal due to numerical error.
849  // Therefore, we use the Approx comparison instead of Equals
850  CHECK_THAT( test, Catch::Approx<double>(truth) );
851  }
852 }
853 
854 
855 TEST_CASE("constant_orthotropic_transient_heat_transfer_material_3d",
856  "[orthotropic],[material],[3D],[heat_transfer],[constant],[transient]")
857 {
858  // Define Material Properties as MAST Parameters
859  MAST::Parameter rho("rho_param", 1234.5); // Density
860  MAST::Parameter cp("cp_param", 908.0); // Specific Heat Capacity
861 
862  // Create field functions to dsitribute these constant parameters throughout the model
863  MAST::ConstantFieldFunction rho_f("rho", rho);
864  MAST::ConstantFieldFunction cp_f("cp", cp);
865 
866  // Initialize the material
868 
869  // Add the material property constant field functions to the material card
870  material.add(rho_f);
871  material.add(cp_f);
872 
873  // Ensure that the material properties were added before doing other checks
874  REQUIRE( material.contains("rho") );
875  REQUIRE( material.contains("cp") );
876 
877  SECTION("material depends on the parameters that it should")
878  {
879  CHECK( material.depends_on(rho) );
880  CHECK( material.depends_on(cp) );
881  }
882 
883  SECTION("material does not depend on other parameters")
884  {
885  MAST::Parameter dummy("dummy", 1.0);
886  CHECK_FALSE( material.depends_on(dummy) );
887  }
888 
889  SECTION("2D thermal capacitance matrix")
890  {
891  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
892  material.capacitance_matrix(3);
893 
894  const libMesh::Point point(2.3, 3.1, 5.2);
895  const Real time = 2.34;
896  RealMatrixX D_cp;
897  mat_capacit(point, time, D_cp);
898 
899  // Hard-coded values for thermal conductivity matrix
900  RealMatrixX D_cp_true = RealMatrixX::Identity(1,1);
901  D_cp_true *= (908.0 * 1234.5);
902 
903  // Convert the test and truth Eigen::Matrix objects to std::vector
904  // since Catch2 has built in methods to compare vectors
905  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_cp);
906  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_cp_true);
907 
908  // Floating point approximations are diffcult to compare since the
909  // values typically aren't exactly equal due to numerical error.
910  // Therefore, we use the Approx comparison instead of Equals
911  CHECK_THAT( test, Catch::Approx<double>(truth) );
912  }
913 }
914 
915 
916 TEST_CASE("constant_orthotropic_dynamic_material_3d",
917  "[orthotropic],[material],[3D],[dynamic],[constant]")
918 {
919  // Define Material Properties as MAST Parameters
920  MAST::Parameter rho("rho_param", 1234.5); // Density
921 
922  // Create field functions to dsitribute these constant parameters throughout the model
923  MAST::ConstantFieldFunction rho_f("rho", rho);
924 
925  // Initialize the material
927 
928  // Add the material property constant field functions to the material card
929  material.add(rho_f);
930 
931  // Ensure that the material properties were added before doing other checks
932  REQUIRE( material.contains("rho") );
933 
934  SECTION("material depends on the parameters that it should")
935  {
936  CHECK( material.depends_on(rho) );
937  }
938 
939  SECTION("material does not depend on other parameters")
940  {
941  MAST::Parameter dummy("dummy", 1.0);
942  CHECK_FALSE( material.depends_on(dummy) );
943  }
944 
945  SECTION("3D inertia matrix")
946  {
947  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
948  material.inertia_matrix(3);
949 
950  const libMesh::Point point(2.3, 3.1, 5.2);
951  const Real time = 2.34;
952  RealMatrixX D_inertia;
953  mat_capacit(point, time, D_inertia);
954 
955  // Hard-coded values for thermal conductivity matrix
956  RealMatrixX D_inertia_true = RealMatrixX::Identity(3,3);
957  D_inertia_true *= 1234.5;
958 
959  // Convert the test and truth Eigen::Matrix objects to std::vector
960  // since Catch2 has built in methods to compare vectors
961  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_inertia);
962  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_inertia_true);
963 
964  // Floating point approximations are diffcult to compare since the
965  // values typically aren't exactly equal due to numerical error.
966  // Therefore, we use the Approx comparison instead of Equals
967  CHECK_THAT( test, Catch::Approx<double>(truth) );
968  }
969 }
TEST_CASE("constant_orthotropic_thermoelastic_material_1d", "[orthotropic],[material],[constant],[1D],[thermoelastic]")
virtual const MAST::FieldFunction< RealMatrixX > & conductance_matrix(const unsigned int dim)
virtual const MAST::FieldFunction< RealMatrixX > & thermal_expansion_matrix(const unsigned int dim)
This is a scalar function whose value can be changed and one that can be used as a design variable in...
Definition: parameter.h:35
virtual const MAST::FieldFunction< RealMatrixX > & inertia_matrix(const unsigned int dim)
libMesh::Real Real
virtual const MAST::FieldFunction< RealMatrixX > & stiffness_matrix(const unsigned int dim, const bool plane_stress=true)
void add(MAST::FunctionBase &f)
adds the function to this card and returns a reference to it.
std::vector< double > eigen_matrix_to_std_vector(RealMatrixX M)
Converts an Eigen Matrix object to a std::vector.
bool contains(const std::string &nm) const
checks if the card contains the specified property value
Matrix< Real, Dynamic, Dynamic > RealMatrixX
Matrix< Real, Dynamic, 1 > RealVectorX
virtual const MAST::FieldFunction< RealMatrixX > & capacitance_matrix(const unsigned int dim)
virtual const MAST::FieldFunction< RealMatrixX > & transverse_shear_stiffness_matrix()
virtual bool depends_on(const MAST::FunctionBase &f) const
returns true if the property card depends on the function f