MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_isotropic_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 
18 TEST_CASE("constant_isotropic_thermoelastic_material_1d",
19  "[isotropic],[material],[constant],[1D],[thermoelastic]")
20 {
21  MAST::Parameter alpha("alpha_param", 5.43e-05); // Coefficient of thermal expansion
22  MAST::ConstantFieldFunction alpha_f("alpha_expansion", alpha);
23 
24  // Initialize the material
26 
27  material.add(alpha_f);
28 
29  REQUIRE( material.contains("alpha_expansion") );
30 
31  SECTION("1D material thermal expansion matrix")
32  {
33  const MAST::FieldFunction<RealMatrixX>& mat_texp =
34  material.thermal_expansion_matrix(1);
35 
36  const libMesh::Point point(2.3, 3.1, 5.2);
37  const Real time = 2.34;
38  RealMatrixX D_texp;
39  mat_texp(point, time, D_texp);
40 
41  // Hard-coded in the true value of material thermal expansion matrix
42  RealMatrixX D_texp_true = RealMatrixX::Zero(2,1);
43  D_texp_true(0,0) = 5.43e-05;
44 
45  // Convert the test and truth Eigen::Matrix objects to std::vector
46  // since Catch2 has built in methods to compare vectors
47  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_texp);
48  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_texp_true);
49 
50  // Floating point approximations are diffcult to compare since the
51  // values typically aren't exactly equal due to numerical error.
52  // Therefore, we use the Approx comparison instead of Equals
53  CHECK_THAT( test, Catch::Approx<double>(truth) );
54  }
55 }
56 
57 
58 TEST_CASE("constant_isotropic_structural_material_1d",
59  "[isotropic],[material],[constant],[1D],[structural]")
60 {
61  // Define Material Properties as MAST Parameters
62  MAST::Parameter E("E_param", 72.0e9); // Modulus of Elasticity
63  MAST::Parameter nu("nu_param", 0.33); // Poisson's ratio
64 
65  // Create field functions to dsitribute these constant parameters throughout the model
66  MAST::ConstantFieldFunction E_f("E", E);
67  MAST::ConstantFieldFunction nu_f("nu", nu);
68 
69  // Initialize the material
71 
72  // Add the material property constant field functions to the material card
73  material.add(E_f);
74  material.add(nu_f);
75 
76  // Ensure that the material properties were added before doing other checks
77  REQUIRE( material.contains("E") );
78  REQUIRE( material.contains("nu") );
79 
80  SECTION("1D material constitutive matrix")
81  {
82  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
83  material.stiffness_matrix(1);
84 
85  const libMesh::Point point(2.3, 3.1, 5.2);
86  const Real time = 2.34;
87  RealMatrixX D;
88  mat_stiffness(point, time, D);
89 
90  // Hard-code in the true value of the material stiffness
91  RealMatrixX D_true = RealMatrixX::Zero(2,2);
92  D_true(0,0) = 72.0e9;
93  D_true(1,1) = 2.706766917293233e+10;
94 
95  // Convert the test and truth Eigen::Matrix objects to std::vector
96  // since Catch2 has built in methods to compare vectors
97  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D);
98  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_true);
99 
100  // Floating point approximations are diffcult to compare since the
101  // values typically aren't exactly equal due to numerical error.
102  // Therefore, we use the Approx comparison instead of Equals
103  CHECK_THAT( test, Catch::Approx<double>(truth) );
104  }
105 }
106 
107 
108 TEST_CASE("constant_isotropic_heat_transfer_material_1d",
109  "[isotropic],[material],[1D],[heat_transfer][constant]")
110 {
111  // Define Material Properties as MAST Parameters
112  MAST::Parameter k("k_param", 237.0); // Thermal Conductivity
113 
114  // Create field functions to dsitribute these constant parameters throughout the model
115  MAST::ConstantFieldFunction k_f("k_th", k);
116 
117  // Initialize the material
119 
120  // Add the material property constant field functions to the material card
121  material.add(k_f);
122 
123  // Ensure that the material properties were added before doing other checks
124  REQUIRE( material.contains("k_th") );
125 
126  SECTION("material depends on the parameters that it should")
127  {
128  CHECK( material.depends_on(k) );
129  }
130 
131  SECTION("material does not depend on other parameters")
132  {
133  MAST::Parameter dummy("dummy", 1.0);
134  CHECK_FALSE( material.depends_on(dummy) );
135  }
136 
137  SECTION("1D thermal conductivity matrix")
138  {
139  const MAST::FieldFunction<RealMatrixX>& mat_conduct =
140  material.conductance_matrix(1);
141 
142  const libMesh::Point point(2.3, 3.1, 5.2);
143  const Real time = 2.34;
144  RealMatrixX D_k;
145  mat_conduct(point, time, D_k);
146 
147  // Hard-coded values for thermal conductivity matrix
148  RealMatrixX D_k_true = RealMatrixX::Identity(1,1);
149  D_k_true *= 237.0;
150 
151  // Convert the test and truth Eigen::Matrix objects to std::vector
152  // since Catch2 has built in methods to compare vectors
153  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_k);
154  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_k_true);
155 
156  // Floating point approximations are diffcult to compare since the
157  // values typically aren't exactly equal due to numerical error.
158  // Therefore, we use the Approx comparison instead of Equals
159  CHECK_THAT( test, Catch::Approx<double>(truth) );
160  }
161 }
162 
163 
164 TEST_CASE("constant_isotropic_transient_heat_transfer_material_1d",
165  "[isotropic],[material],[1D],[heat_transfer],[constant],[transient]")
166 {
167  // Define Material Properties as MAST Parameters
168  MAST::Parameter rho("rho_param", 1234.5); // Density
169  MAST::Parameter cp("cp_param", 908.0); // Specific Heat Capacity
170 
171  // Create field functions to dsitribute these constant parameters throughout the model
172  MAST::ConstantFieldFunction rho_f("rho", rho);
173  MAST::ConstantFieldFunction cp_f("cp", cp);
174 
175  // Initialize the material
177 
178  // Add the material property constant field functions to the material card
179  material.add(rho_f);
180  material.add(cp_f);
181 
182  // Ensure that the material properties were added before doing other checks
183  REQUIRE( material.contains("rho") );
184  REQUIRE( material.contains("cp") );
185 
186  SECTION("material depends on the parameters that it should")
187  {
188  CHECK( material.depends_on(rho) );
189  CHECK( material.depends_on(cp) );
190  }
191 
192  SECTION("material does not depend on other parameters")
193  {
194  MAST::Parameter dummy("dummy", 1.0);
195  CHECK_FALSE( material.depends_on(dummy) );
196  }
197 
198  SECTION("1D thermal capacitance matrix")
199  {
200  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
201  material.capacitance_matrix(1);
202 
203  const libMesh::Point point(2.3, 3.1, 5.2);
204  const Real time = 2.34;
205  RealMatrixX D_cp;
206  mat_capacit(point, time, D_cp);
207 
208  // Hard-coded values for thermal conductivity matrix
209  RealMatrixX D_cp_true = RealMatrixX::Identity(1,1);
210  D_cp_true *= (908.0 * 1234.5);
211 
212  // Convert the test and truth Eigen::Matrix objects to std::vector
213  // since Catch2 has built in methods to compare vectors
214  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_cp);
215  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_cp_true);
216 
217  // Floating point approximations are diffcult to compare since the
218  // values typically aren't exactly equal due to numerical error.
219  // Therefore, we use the Approx comparison instead of Equals
220  CHECK_THAT( test, Catch::Approx<double>(truth) );
221  }
222 }
223 
224 
225 TEST_CASE("constant_isotropic_thermoelastic_material_2d",
226  "[isotropic],[material],[2D],[thermoelastic][constant]")
227 {
228  MAST::Parameter alpha("alpha_param", 5.43e-05); // Coefficient of thermal expansion
229  MAST::ConstantFieldFunction alpha_f("alpha_expansion", alpha);
230 
231  // Initialize the material
233 
234  material.add(alpha_f);
235 
236  REQUIRE( material.contains("alpha_expansion") );
237 
238  SECTION("material does not depend on other parameters")
239  {
240  MAST::Parameter dummy("dummy", 1.0);
241  CHECK_FALSE( material.depends_on(dummy) );
242  }
243 
244  SECTION("2D material thermal expansion matrix")
245  {
246  const MAST::FieldFunction<RealMatrixX>& mat_texp =
247  material.thermal_expansion_matrix(2);
248 
249  const libMesh::Point point(2.3, 3.1, 5.2);
250  const Real time = 2.34;
251  RealMatrixX D_texp;
252  mat_texp(point, time, D_texp);
253 
254  // Hard-coded in the true value of material thermal expansion matrix
255  RealMatrixX D_texp_true = RealMatrixX::Zero(3,1);
256  D_texp_true(0,0) = 5.43e-05;
257  D_texp_true(1,0) = 5.43e-05;
258 
259  // Convert the test and truth Eigen::Matrix objects to std::vector
260  // since Catch2 has built in methods to compare vectors
261  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_texp);
262  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_texp_true);
263 
264  // Floating point approximations are diffcult to compare since the
265  // values typically aren't exactly equal due to numerical error.
266  // Therefore, we use the Approx comparison instead of Equals
267  CHECK_THAT( test, Catch::Approx<double>(truth) );
268  }
269 }
270 
271 
272 TEST_CASE("constant_isotropic_structural_material_2d",
273  "[isotropic],[material],[2D],[structural],[constant]")
274 {
275  // Define Material Properties as MAST Parameters
276  MAST::Parameter E("E_param", 72.0e9); // Modulus of Elasticity
277  MAST::Parameter nu("nu_param", 0.33); // Poisson's ratio
278 
279  const Real G = E() / (2.0 * (1.0+nu()));
280 
281  // Create field functions to dsitribute these constant parameters throughout the model
282  MAST::ConstantFieldFunction E_f("E", E);
283  MAST::ConstantFieldFunction nu_f("nu", nu);
284 
285  // Initialize the material
287 
288  // Add the material property constant field functions to the material card
289  material.add(E_f);
290  material.add(nu_f);
291 
292  // Ensure that the material properties were added before doing other checks
293  REQUIRE( material.contains("E") );
294  REQUIRE( material.contains("nu") );
295 
296  SECTION("material depends on the parameters that it should")
297  {
298  CHECK( material.depends_on(E) );
299  CHECK( material.depends_on(nu) );
300  }
301 
302  SECTION("material does not depend on other parameters")
303  {
304  MAST::Parameter dummy("dummy", 1.0);
305  CHECK_FALSE( material.depends_on(dummy) );
306  }
307 
308  SECTION("2D plane stress material constitutive matrix is correct")
309  {
310  const bool is_plane_stress = true;
311  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
312  material.stiffness_matrix(2, is_plane_stress);
313 
314  const libMesh::Point point(2.3, 3.1, 5.2);
315  const Real time = 2.34;
316  RealMatrixX K_mat;
317  mat_stiffness(point, time, K_mat);
318 
319  // Check for symmetry
320 
321  // Convert the test and truth Eigen::Matrix objects to std::vector
322  // since Catch2 has built in methods to compare vectors
323  std::vector<double> test = TEST::eigen_matrix_to_std_vector(K_mat.transpose());
324  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(K_mat);
325 
326  // Floating point approximations are diffcult to compare since the
327  // values typically aren't exactly equal due to numerical error.
328  // Therefore, we use the Approx comparison instead of Equals
329  REQUIRE_THAT( test, Catch::Approx<double>(truth) );
330 
331  // Check for positive definiteness
332  SelfAdjointEigenSolver<RealMatrixX> eigensolver(K_mat);
333  if (eigensolver.info() != Success)
334  {
335  std::cout << "eigensolver failed to converge!" << std::endl;
336  REQUIRE(false); // Force test failure.
337  }
338  RealVectorX eigenvalues = eigensolver.eigenvalues();
339  for (uint i=0; i<K_mat.cols(); i++)
340  {
341  REQUIRE( eigenvalues(i) > 0.0 );
342  }
343 
344  // Hard-code in the true value of the material stiffness
345  RealMatrixX K_mat_true = RealMatrixX::Zero(3,3);
346  K_mat_true(0,0) = 8.079901245651442e+10;
347  K_mat_true(1,1) = 8.079901245651442e+10;
348  K_mat_true(0,1) = 2.666367411064976e+10;
349  K_mat_true(1,0) = 2.666367411064976e+10;
350  K_mat_true(2,2) = 2.706766917293233e+10;
351 
352  // Convert the test and truth Eigen::Matrix objects to std::vector
353  // since Catch2 has built in methods to compare vectors
354  test = TEST::eigen_matrix_to_std_vector(K_mat);
355  truth = TEST::eigen_matrix_to_std_vector(K_mat_true);
356 
357  // Floating point approximations are diffcult to compare since the
358  // values typically aren't exactly equal due to numerical error.
359  // Therefore, we use the Approx comparison instead of Equals
360  CHECK_THAT( test, Catch::Approx<double>(truth) );
361  }
362 
363  SECTION("2D plane strain material constitutive matrix is correct")
364  {
365  const bool is_plane_stress = false;
366  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
367  material.stiffness_matrix(2, is_plane_stress);
368 
369  const libMesh::Point point(2.3, 3.1, 5.2);
370  const Real time = 2.34;
371  RealMatrixX K_mat;
372  mat_stiffness(point, time, K_mat);
373 
374  // Check for symmetry
375 
376  // Convert the test and truth Eigen::Matrix objects to std::vector
377  // since Catch2 has built in methods to compare vectors
378  std::vector<double> test = TEST::eigen_matrix_to_std_vector(K_mat.transpose());
379  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(K_mat);
380 
381  // Floating point approximations are diffcult to compare since the
382  // values typically aren't exactly equal due to numerical error.
383  // Therefore, we use the Approx comparison instead of Equals
384  REQUIRE_THAT( test, Catch::Approx<double>(truth) );
385 
386  // Check for positive definiteness
387  SelfAdjointEigenSolver<RealMatrixX> eigensolver(K_mat);
388  if (eigensolver.info() != Success)
389  {
390  std::cout << "eigensolver failed to converge!" << std::endl;
391  REQUIRE(false); // Force test failure.
392  }
393  RealVectorX eigenvalues = eigensolver.eigenvalues();
394  for (uint i=0; i<K_mat.cols(); i++)
395  {
396  REQUIRE( eigenvalues(i) > 0.0 );
397  }
398 
399  // Hard-code in the true value of the material stiffness
400  RealMatrixX K_mat_true = RealMatrixX::Zero(3,3);
401  K_mat_true(0,0) = 1.066784608580274e+11;
402  K_mat_true(1,1) = 1.066784608580274e+11;
403  K_mat_true(0,1) = 0.525431225121628e+11;
404  K_mat_true(1,0) = 0.525431225121628e+11;
405  K_mat_true(2,2) = 0.270676691729323e+11;
406 
407  // Convert the test and truth Eigen::Matrix objects to std::vector
408  // since Catch2 has built in methods to compare vectors
409  test = TEST::eigen_matrix_to_std_vector(K_mat);
410  truth = TEST::eigen_matrix_to_std_vector(K_mat_true);
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  CHECK_THAT( test, Catch::Approx<double>(truth) );
416  }
417 
418  SECTION("transverse shear material stiffness matrix")
419  {
420  const MAST::FieldFunction<RealMatrixX>& mat_trans_shear =
422 
423  const libMesh::Point point(2.3, 3.1, 5.2);
424  const Real time = 2.34;
425  RealMatrixX D_trans_shear;
426  mat_trans_shear(point, time, D_trans_shear);
427 
428  // Hard-coded in the true value of the material transverse shear stiffness
429  RealMatrixX D_trans_shear_true = RealMatrixX::Zero(2,2);
430  D_trans_shear_true(0,0) = G;//2.255639097744361e+10;
431  D_trans_shear_true(1,1) = G;//2.255639097744361e+10;
432 
433  // Convert the test and truth Eigen::Matrix objects to std::vector
434  // since Catch2 has built in methods to compare vectors
435  std::vector<double> test(D_trans_shear.data(), D_trans_shear.data()+D_trans_shear.rows()*D_trans_shear.cols());
436  std::vector<double> truth(D_trans_shear_true.data(), D_trans_shear_true.data()+D_trans_shear_true.rows()*D_trans_shear_true.cols());
437 
438  // Floating point approximations are diffcult to compare since the
439  // values typically aren't exactly equal due to numerical error.
440  // Therefore, we use the Approx comparison instead of Equals
441  CHECK_THAT( test, Catch::Approx<double>(truth) );
442  }
443 }
444 
445 
446 TEST_CASE("constant_isotropic_heat_transfer_material_2d",
447  "[isotropic],[material],[2D],[heat_transfer][constant]")
448 {
449  // Define Material Properties as MAST Parameters
450  MAST::Parameter k("k_param", 237.0); // Thermal Conductivity
451 
452  // Create field functions to dsitribute these constant parameters throughout the model
453  MAST::ConstantFieldFunction k_f("k_th", k);
454 
455  // Initialize the material
457 
458  // Add the material property constant field functions to the material card
459  material.add(k_f);
460 
461  // Ensure that the material properties were added before doing other checks
462  REQUIRE( material.contains("k_th") );
463 
464  SECTION("material depends on the parameters that it should")
465  {
466  CHECK( material.depends_on(k) );
467  }
468 
469  SECTION("material does not depend on other parameters")
470  {
471  MAST::Parameter dummy("dummy", 1.0);
472  CHECK_FALSE( material.depends_on(dummy) );
473  }
474 
475  SECTION("2D thermal conductivity matrix")
476  {
477  const MAST::FieldFunction<RealMatrixX>& mat_conduct =
478  material.conductance_matrix(2);
479 
480  const libMesh::Point point(2.3, 3.1, 5.2);
481  const Real time = 2.34;
482  RealMatrixX D_k;
483  mat_conduct(point, time, D_k);
484 
485  // Hard-coded values for thermal conductivity matrix
486  RealMatrixX D_k_true = RealMatrixX::Identity(2,2);
487  D_k_true *= 237.0;
488 
489  // Convert the test and truth Eigen::Matrix objects to std::vector
490  // since Catch2 has built in methods to compare vectors
491  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_k);
492  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_k_true);
493 
494  // Floating point approximations are diffcult to compare since the
495  // values typically aren't exactly equal due to numerical error.
496  // Therefore, we use the Approx comparison instead of Equals
497  CHECK_THAT( test, Catch::Approx<double>(truth) );
498  }
499 }
500 
501 
502 TEST_CASE("constant_isotropic_transient_heat_transfer_material_2d",
503  "[isotropic],[material],[2D],[heat_transfer],[constant],[transient]")
504 {
505  // Define Material Properties as MAST Parameters
506  MAST::Parameter rho("rho_param", 1234.5); // Density
507  MAST::Parameter cp("cp_param", 908.0); // Specific Heat Capacity
508 
509  // Create field functions to dsitribute these constant parameters throughout the model
510  MAST::ConstantFieldFunction rho_f("rho", rho);
511  MAST::ConstantFieldFunction cp_f("cp", cp);
512 
513  // Initialize the material
515 
516  // Add the material property constant field functions to the material card
517  material.add(rho_f);
518  material.add(cp_f);
519 
520  // Ensure that the material properties were added before doing other checks
521  REQUIRE( material.contains("rho") );
522  REQUIRE( material.contains("cp") );
523 
524  SECTION("material depends on the parameters that it should")
525  {
526  CHECK( material.depends_on(rho) );
527  CHECK( material.depends_on(cp) );
528  }
529 
530  SECTION("material does not depend on other parameters")
531  {
532  MAST::Parameter dummy("dummy", 1.0);
533  CHECK_FALSE( material.depends_on(dummy) );
534  }
535 
536  SECTION("2D thermal capacitance matrix")
537  {
538  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
539  material.capacitance_matrix(2);
540 
541  const libMesh::Point point(2.3, 3.1, 5.2);
542  const Real time = 2.34;
543  RealMatrixX D_cp;
544  mat_capacit(point, time, D_cp);
545 
546  // Hard-coded values for thermal conductivity matrix
547  RealMatrixX D_cp_true = RealMatrixX::Identity(1,1);
548  D_cp_true *= (908.0 * 1234.5);
549 
550  // Convert the test and truth Eigen::Matrix objects to std::vector
551  // since Catch2 has built in methods to compare vectors
552  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_cp);
553  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_cp_true);
554 
555  // Floating point approximations are diffcult to compare since the
556  // values typically aren't exactly equal due to numerical error.
557  // Therefore, we use the Approx comparison instead of Equals
558  CHECK_THAT( test, Catch::Approx<double>(truth) );
559  }
560 }
561 
562 
563 TEST_CASE("constant_isotropic_thermoelastic_material_3d",
564  "[isotropic],[material],[3D],[thermoelastic][constant]")
565 {
566  MAST::Parameter alpha("alpha_param", 5.43e-05); // Coefficient of thermal expansion
567  MAST::ConstantFieldFunction alpha_f("alpha_expansion", alpha);
568 
569  // Initialize the material
571 
572  material.add(alpha_f);
573 
574  REQUIRE( material.contains("alpha_expansion") );
575 
576  SECTION("material does not depend on other parameters")
577  {
578  MAST::Parameter dummy("dummy", 1.0);
579  CHECK_FALSE( material.depends_on(dummy) );
580  }
581 
582  SECTION("3D material thermal expansion matrix")
583  {
584  const MAST::FieldFunction<RealMatrixX>& mat_texp =
585  material.thermal_expansion_matrix(3);
586 
587  const libMesh::Point point(2.3, 3.1, 5.2);
588  const Real time = 2.34;
589  RealMatrixX D_texp;
590  mat_texp(point, time, D_texp);
591 
592  // Hard-coded in the true value of material thermal expansion matrix
593  RealMatrixX D_texp_true = RealMatrixX::Zero(6,1);
594  D_texp_true(0,0) = 5.43e-05;
595  D_texp_true(1,0) = 5.43e-05;
596  D_texp_true(2,0) = 5.43e-05;
597 
598  // Convert the test and truth Eigen::Matrix objects to std::vector
599  // since Catch2 has built in methods to compare vectors
600  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_texp);
601  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_texp_true);
602 
603  // Floating point approximations are diffcult to compare since the
604  // values typically aren't exactly equal due to numerical error.
605  // Therefore, we use the Approx comparison instead of Equals
606  CHECK_THAT( test, Catch::Approx<double>(truth) );
607  }
608 }
609 
610 
611 TEST_CASE("constant_isotropic_structural_material_3d",
612  "[isotropic],[material],[3D],[structural],[constant]")
613 {
614  // Define Material Properties as MAST Parameters
615  MAST::Parameter E("E_param", 72.0e9); // Modulus of Elasticity
616  MAST::Parameter nu("nu_param", 0.33); // Poisson's ratio
617 
618  const Real G = E() / (2.0 * (1.0+nu()));
619 
620  // Create field functions to dsitribute these constant parameters throughout the model
621  MAST::ConstantFieldFunction E_f("E", E);
622  MAST::ConstantFieldFunction nu_f("nu", nu);
623 
624  // Initialize the material
626 
627  // Add the material property constant field functions to the material card
628  material.add(E_f);
629  material.add(nu_f);
630 
631  // Ensure that the material properties were added before doing other checks
632  REQUIRE( material.contains("E") );
633  REQUIRE( material.contains("nu") );
634 
635  SECTION("material depends on the parameters that it should")
636  {
637  CHECK( material.depends_on(E) );
638  CHECK( material.depends_on(nu) );
639  }
640 
641  SECTION("material does not depend on other parameters")
642  {
643  MAST::Parameter dummy("dummy", 1.0);
644  CHECK_FALSE( material.depends_on(dummy) );
645  }
646 
647  SECTION("3D material constitutive matrix is correct")
648  {
649  const bool is_plane_stress = true;
650  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
651  material.stiffness_matrix(3);
652 
653  const libMesh::Point point(2.3, 3.1, 5.2);
654  const Real time = 2.34;
655  RealMatrixX K_mat;
656  mat_stiffness(point, time, K_mat);
657 
658  // Check for symmetry
659 
660  // Convert the test and truth Eigen::Matrix objects to std::vector
661  // since Catch2 has built in methods to compare vectors
662  std::vector<double> test = TEST::eigen_matrix_to_std_vector(K_mat.transpose());
663  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(K_mat);
664 
665  // Floating point approximations are diffcult to compare since the
666  // values typically aren't exactly equal due to numerical error.
667  // Therefore, we use the Approx comparison instead of Equals
668  REQUIRE_THAT( test, Catch::Approx<double>(truth) );
669 
670  // Check for positive definiteness
671  SelfAdjointEigenSolver<RealMatrixX> eigensolver(K_mat);
672  if (eigensolver.info() != Success)
673  {
674  std::cout << "eigensolver failed to converge!" << std::endl;
675  REQUIRE(false); // Force test failure.
676  }
677  RealVectorX eigenvalues = eigensolver.eigenvalues();
678  for (uint i=0; i<K_mat.cols(); i++)
679  {
680  REQUIRE( eigenvalues(i) > 0.0 );
681  }
682 
683  // Hard-code in the true value of the material stiffness
684  RealMatrixX K_mat_true = RealMatrixX::Zero(6,6);
685  K_mat_true(0,0) = 1.066784608580274e+11;
686  K_mat_true(1,1) = 1.066784608580274e+11;
687  K_mat_true(2,2) = 1.066784608580274e+11;
688  K_mat_true(3,3) = 0.541353383458647e+11/2.;
689  K_mat_true(4,4) = 0.541353383458647e+11/2.;
690  K_mat_true(5,5) = 0.541353383458647e+11/2.;
691 
692  K_mat_true(0,1) = K_mat_true(1,0) = 0.525431225121628e+11;
693  K_mat_true(0,2) = K_mat_true(2,0) = 0.525431225121628e+11;
694  K_mat_true(1,2) = K_mat_true(2,1) = 0.525431225121628e+11;
695 
696 
697  // Convert the test and truth Eigen::Matrix objects to std::vector
698  // since Catch2 has built in methods to compare vectors
699  test = TEST::eigen_matrix_to_std_vector(K_mat);
700  truth = TEST::eigen_matrix_to_std_vector(K_mat_true);
701 
702  // Floating point approximations are diffcult to compare since the
703  // values typically aren't exactly equal due to numerical error.
704  // Therefore, we use the Approx comparison instead of Equals
705  CHECK_THAT( test, Catch::Approx<double>(truth) );
706  }
707 }
708 
709 
710 TEST_CASE("constant_isotropic_heat_transfer_material_3d",
711  "[isotropic],[material],[3D],[heat_transfer][constant]")
712 {
713  // Define Material Properties as MAST Parameters
714  MAST::Parameter k("k_param", 237.0); // Thermal Conductivity
715 
716  // Create field functions to dsitribute these constant parameters throughout the model
717  MAST::ConstantFieldFunction k_f("k_th", k);
718 
719  // Initialize the material
721 
722  // Add the material property constant field functions to the material card
723  material.add(k_f);
724 
725  // Ensure that the material properties were added before doing other checks
726  REQUIRE( material.contains("k_th") );
727 
728  SECTION("material depends on the parameters that it should")
729  {
730  CHECK( material.depends_on(k) );
731  }
732 
733  SECTION("material does not depend on other parameters")
734  {
735  MAST::Parameter dummy("dummy", 1.0);
736  CHECK_FALSE( material.depends_on(dummy) );
737  }
738 
739  SECTION("2D thermal conductivity matrix")
740  {
741  const MAST::FieldFunction<RealMatrixX>& mat_conduct =
742  material.conductance_matrix(3);
743 
744  const libMesh::Point point(2.3, 3.1, 5.2);
745  const Real time = 2.34;
746  RealMatrixX D_k;
747  mat_conduct(point, time, D_k);
748 
749  // Hard-coded values for thermal conductivity matrix
750  RealMatrixX D_k_true = RealMatrixX::Identity(3,3);
751  D_k_true *= 237.0;
752 
753  // Convert the test and truth Eigen::Matrix objects to std::vector
754  // since Catch2 has built in methods to compare vectors
755  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_k);
756  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_k_true);
757 
758  // Floating point approximations are diffcult to compare since the
759  // values typically aren't exactly equal due to numerical error.
760  // Therefore, we use the Approx comparison instead of Equals
761  CHECK_THAT( test, Catch::Approx<double>(truth) );
762  }
763 }
764 
765 
766 TEST_CASE("constant_isotropic_transient_heat_transfer_material_3d",
767  "[isotropic],[material],[3D],[heat_transfer],[constant],[transient]")
768 {
769  // Define Material Properties as MAST Parameters
770  MAST::Parameter rho("rho_param", 1234.5); // Density
771  MAST::Parameter cp("cp_param", 908.0); // Specific Heat Capacity
772 
773  // Create field functions to dsitribute these constant parameters throughout the model
774  MAST::ConstantFieldFunction rho_f("rho", rho);
775  MAST::ConstantFieldFunction cp_f("cp", cp);
776 
777  // Initialize the material
779 
780  // Add the material property constant field functions to the material card
781  material.add(rho_f);
782  material.add(cp_f);
783 
784  // Ensure that the material properties were added before doing other checks
785  REQUIRE( material.contains("rho") );
786  REQUIRE( material.contains("cp") );
787 
788  SECTION("material depends on the parameters that it should")
789  {
790  CHECK( material.depends_on(rho) );
791  CHECK( material.depends_on(cp) );
792  }
793 
794  SECTION("material does not depend on other parameters")
795  {
796  MAST::Parameter dummy("dummy", 1.0);
797  CHECK_FALSE( material.depends_on(dummy) );
798  }
799 
800  SECTION("2D thermal capacitance matrix")
801  {
802  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
803  material.capacitance_matrix(3);
804 
805  const libMesh::Point point(2.3, 3.1, 5.2);
806  const Real time = 2.34;
807  RealMatrixX D_cp;
808  mat_capacit(point, time, D_cp);
809 
810  // Hard-coded values for thermal conductivity matrix
811  RealMatrixX D_cp_true = RealMatrixX::Identity(1,1);
812  D_cp_true *= (908.0 * 1234.5);
813 
814  // Convert the test and truth Eigen::Matrix objects to std::vector
815  // since Catch2 has built in methods to compare vectors
816  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_cp);
817  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_cp_true);
818 
819  // Floating point approximations are diffcult to compare since the
820  // values typically aren't exactly equal due to numerical error.
821  // Therefore, we use the Approx comparison instead of Equals
822  CHECK_THAT( test, Catch::Approx<double>(truth) );
823  }
824 }
825 
826 
827 TEST_CASE("constant_isotropic_dynamic_material_3d",
828  "[isotropic],[material],[3D],[dynamic],[constant]")
829 {
830  // Define Material Properties as MAST Parameters
831  MAST::Parameter rho("rho_param", 1234.5); // Density
832 
833  // Create field functions to dsitribute these constant parameters throughout the model
834  MAST::ConstantFieldFunction rho_f("rho", rho);
835 
836  // Initialize the material
838 
839  // Add the material property constant field functions to the material card
840  material.add(rho_f);
841 
842  // Ensure that the material properties were added before doing other checks
843  REQUIRE( material.contains("rho") );
844 
845  SECTION("material depends on the parameters that it should")
846  {
847  CHECK( material.depends_on(rho) );
848  }
849 
850  SECTION("material does not depend on other parameters")
851  {
852  MAST::Parameter dummy("dummy", 1.0);
853  CHECK_FALSE( material.depends_on(dummy) );
854  }
855 
856  SECTION("3D inertia matrix")
857  {
858  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
859  material.inertia_matrix(3);
860 
861  const libMesh::Point point(2.3, 3.1, 5.2);
862  const Real time = 2.34;
863  RealMatrixX D_inertia;
864  mat_capacit(point, time, D_inertia);
865 
866  // Hard-coded values for thermal conductivity matrix
867  RealMatrixX D_inertia_true = RealMatrixX::Identity(3,3);
868  D_inertia_true *= 1234.5;
869 
870  // Convert the test and truth Eigen::Matrix objects to std::vector
871  // since Catch2 has built in methods to compare vectors
872  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_inertia);
873  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_inertia_true);
874 
875  // Floating point approximations are diffcult to compare since the
876  // values typically aren't exactly equal due to numerical error.
877  // Therefore, we use the Approx comparison instead of Equals
878  CHECK_THAT( test, Catch::Approx<double>(truth) );
879  }
880 }
virtual const MAST::FieldFunction< RealMatrixX > & inertia_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
libMesh::Real Real
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
virtual const MAST::FieldFunction< RealMatrixX > & conductance_matrix(const unsigned int dim)
Matrix< Real, Dynamic, 1 > RealVectorX
virtual const MAST::FieldFunction< RealMatrixX > & capacitance_matrix(const unsigned int dim)
virtual const MAST::FieldFunction< RealMatrixX > & stiffness_matrix(const unsigned int dim, const bool plane_stress=true)
TEST_CASE("constant_isotropic_thermoelastic_material_1d", "[isotropic],[material],[constant],[1D],[thermoelastic]")
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