MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_isotropic_material_sensitivity.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 // TODO: Check this. May need to remove the effect of Kappa from the calculations.
21 
22 TEST_CASE("constant_isotropic_thermoelastic_material_1d_sensitivity",
23  "[isotropic],[material],[constant],[1D],[thermoelastic],[sensitivity]")
24 {
25  MAST::Parameter alpha("alpha_param", 5.43e-05); // Coefficient of thermal expansion
26  MAST::ConstantFieldFunction alpha_f("alpha_expansion", alpha);
27 
28  // Initialize the material
30 
31  material.add(alpha_f);
32 
33  REQUIRE( material.contains("alpha_expansion") );
34 
35  SECTION("1D material thermal expansion matrix derivative w.r.t. coefficient of thermal expansion")
36  {
37  const MAST::FieldFunction<RealMatrixX>& mat_texp =
38  material.thermal_expansion_matrix(1);
39 
40  const libMesh::Point point(2.3, 3.1, 5.2);
41  const Real time = 2.34;
42  RealMatrixX dD_texp;
43  mat_texp.derivative(alpha, point, time, dD_texp);
44 
45  // Hard-coded in the true value of material thermal expansion matrix
46  RealMatrixX dD_texp_true = RealMatrixX::Zero(2,1);
47  dD_texp_true(0,0) = 1.0;
48 
49  // Convert the test and truth Eigen::Matrix objects to std::vector
50  // since Catch2 has built in methods to compare vectors
51  std::vector<double> test = TEST::eigen_matrix_to_std_vector(dD_texp);
52  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(dD_texp_true);
53 
54  // Floating point approximations are diffcult to compare since the
55  // values typically aren't exactly equal due to numerical error.
56  // Therefore, we use the Approx comparison instead of Equals
57  CHECK_THAT( test, Catch::Approx<double>(truth) );
58  }
59 }
60 
61 
62 TEST_CASE("constant_isotropic_structural_material_1d_sensitivity",
63  "[isotropic],[material],[constant],[1D],[structural],[sensitivity]")
64 {
65  // Define Material Properties as MAST Parameters
66  MAST::Parameter E("E_param", 72.0e9); // Modulus of Elasticity
67  MAST::Parameter nu("nu_param", 0.33); // Poisson's ratio
68  MAST::Parameter kappa("kappa_param", 5.0/6.0); // Shear coefficient
69 
70  // Create field functions to dsitribute these constant parameters throughout the model
71  MAST::ConstantFieldFunction E_f("E", E);
72  MAST::ConstantFieldFunction nu_f("nu", nu);
73  MAST::ConstantFieldFunction kappa_f("kappa", kappa);
74 
75  // Initialize the material
77 
78  // Add the material property constant field functions to the material card
79  material.add(E_f);
80  material.add(nu_f);
81  material.add(kappa_f);
82 
83  // Ensure that the material properties were added before doing other checks
84  REQUIRE( material.contains("E") );
85  REQUIRE( material.contains("nu") );
86  REQUIRE( material.contains("kappa") );
87 
88  SECTION("1D material constitutive matrix derivative w.r.t. elastic modulus")
89  {
90  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
91  material.stiffness_matrix(1);
92 
93  const libMesh::Point point(2.3, 3.1, 5.2);
94  const Real time = 2.34;
95  RealMatrixX dD;
96  mat_stiffness.derivative(E, point, time, dD);
97 
98  // Hard-code in the true value of the material stiffness
99  RealMatrixX dD_true = RealMatrixX::Zero(2,2);
100  dD_true(0,0) = 1.0;
101  dD_true(1,1) = 0.375939849624060;
102 
103  // Convert the test and truth Eigen::Matrix objects to std::vector
104  // since Catch2 has built in methods to compare vectors
105  std::vector<double> test = TEST::eigen_matrix_to_std_vector(dD);
106  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(dD_true);
107 
108  // Floating point approximations are diffcult to compare since the
109  // values typically aren't exactly equal due to numerical error.
110  // Therefore, we use the Approx comparison instead of Equals
111  CHECK_THAT( test, Catch::Approx<double>(truth) );
112  }
113 
114  SECTION("1D material constitutive matrix derivative w.r.t. poisson's ratio")
115  {
116  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
117  material.stiffness_matrix(1);
118 
119  const libMesh::Point point(2.3, 3.1, 5.2);
120  const Real time = 2.34;
121  RealMatrixX dD;
122  mat_stiffness.derivative(nu, point, time, dD);
123 
124  // Hard-code in the true value of the material stiffness
125  RealMatrixX dD_true = RealMatrixX::Zero(2,2);
126  dD_true(0,0) = 0.0;
127  dD_true(1,1) = -0.203516309570920e+11;
128 
129  // Convert the test and truth Eigen::Matrix objects to std::vector
130  // since Catch2 has built in methods to compare vectors
131  std::vector<double> test = TEST::eigen_matrix_to_std_vector(dD);
132  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(dD_true);
133 
134  // Floating point approximations are diffcult to compare since the
135  // values typically aren't exactly equal due to numerical error.
136  // Therefore, we use the Approx comparison instead of Equals
137  CHECK_THAT( test, Catch::Approx<double>(truth) );
138  }
139 }
140 
141 
142 TEST_CASE("constant_isotropic_heat_transfer_material_1d_sensitivity",
143  "[isotropic],[material],[1D],[heat_transfer],[constant],[sensitivity]")
144 {
145  // Define Material Properties as MAST Parameters
146  MAST::Parameter k("k_param", 237.0); // Thermal Conductivity
147 
148  // Create field functions to dsitribute these constant parameters throughout the model
149  MAST::ConstantFieldFunction k_f("k_th", k);
150 
151  // Initialize the material
153 
154  // Add the material property constant field functions to the material card
155  material.add(k_f);
156 
157  // Ensure that the material properties were added before doing other checks
158  REQUIRE( material.contains("k_th") );
159 
160  SECTION("material depends on the parameters that it should")
161  {
162  CHECK( material.depends_on(k) );
163  }
164 
165  SECTION("material does not depend on other parameters")
166  {
167  MAST::Parameter dummy("dummy", 1.0);
168  CHECK_FALSE( material.depends_on(dummy) );
169  }
170 
171  SECTION("1D thermal conductivity matrix derivative w.r.t. thermal conductivity")
172  {
173  const MAST::FieldFunction<RealMatrixX>& mat_conduct =
174  material.conductance_matrix(1);
175 
176  const libMesh::Point point(2.3, 3.1, 5.2);
177  const Real time = 2.34;
178  RealMatrixX dD_k;
179  mat_conduct.derivative(k, point, time, dD_k);
180 
181  // Hard-coded values for thermal conductivity matrix
182  RealMatrixX dD_k_true = RealMatrixX::Identity(1,1);
183 
184  // Convert the test and truth Eigen::Matrix objects to std::vector
185  // since Catch2 has built in methods to compare vectors
186  std::vector<double> test = TEST::eigen_matrix_to_std_vector(dD_k);
187  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(dD_k_true);
188 
189  // Floating point approximations are diffcult to compare since the
190  // values typically aren't exactly equal due to numerical error.
191  // Therefore, we use the Approx comparison instead of Equals
192  CHECK_THAT( test, Catch::Approx<double>(truth) );
193  }
194 }
195 
196 
197 TEST_CASE("constant_isotropic_transient_heat_transfer_material_1d_sensitivity",
198  "[isotropic],[material],[1D],[heat_transfer],[constant],[transient],[sensitivity]")
199 {
200  // Define Material Properties as MAST Parameters
201  MAST::Parameter rho("rho_param", 1234.5); // Density
202  MAST::Parameter cp("cp_param", 908.0); // Specific Heat Capacity
203 
204  // Create field functions to dsitribute these constant parameters throughout the model
205  MAST::ConstantFieldFunction rho_f("rho", rho);
206  MAST::ConstantFieldFunction cp_f("cp", cp);
207 
208  // Initialize the material
210 
211  // Add the material property constant field functions to the material card
212  material.add(rho_f);
213  material.add(cp_f);
214 
215  // Ensure that the material properties were added before doing other checks
216  REQUIRE( material.contains("rho") );
217  REQUIRE( material.contains("cp") );
218 
219  SECTION("material depends on the parameters that it should")
220  {
221  CHECK( material.depends_on(rho) );
222  CHECK( material.depends_on(cp) );
223  }
224 
225  SECTION("material does not depend on other parameters")
226  {
227  MAST::Parameter dummy("dummy", 1.0);
228  CHECK_FALSE( material.depends_on(dummy) );
229  }
230 
231  SECTION("1D thermal capacitance matrix derivative w.r.t. specifc heat capacity")
232  {
233  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
234  material.capacitance_matrix(1);
235 
236  const libMesh::Point point(2.3, 3.1, 5.2);
237  const Real time = 2.34;
238  RealMatrixX dD_cp;
239  mat_capacit.derivative(cp, point, time, dD_cp);
240 
241  // Hard-coded values for thermal conductivity matrix
242  RealMatrixX dD_cp_true = RealMatrixX::Identity(1,1);
243  dD_cp_true *= 1234.5;
244 
245  // Convert the test and truth Eigen::Matrix objects to std::vector
246  // since Catch2 has built in methods to compare vectors
247  std::vector<double> test = TEST::eigen_matrix_to_std_vector(dD_cp);
248  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(dD_cp_true);
249 
250  // Floating point approximations are diffcult to compare since the
251  // values typically aren't exactly equal due to numerical error.
252  // Therefore, we use the Approx comparison instead of Equals
253  CHECK_THAT( test, Catch::Approx<double>(truth) );
254  }
255 }
256 
257 
258 TEST_CASE("constant_isotropic_thermoelastic_material_2d_sensitivity",
259  "[isotropic],[material],[2D],[thermoelastic],[constant],[sensitivity]")
260 {
261  MAST::Parameter alpha("alpha_param", 5.43e-05); // Coefficient of thermal expansion
262  MAST::ConstantFieldFunction alpha_f("alpha_expansion", alpha);
263 
264  // Initialize the material
266 
267  material.add(alpha_f);
268 
269  REQUIRE( material.contains("alpha_expansion") );
270 
271  SECTION("material does not depend on other parameters")
272  {
273  MAST::Parameter dummy("dummy", 1.0);
274  CHECK_FALSE( material.depends_on(dummy) );
275  }
276 
277  SECTION("2D material thermal expansion matrix derivative w.r.t. coefficient of thermal expansion")
278  {
279  const MAST::FieldFunction<RealMatrixX>& mat_texp =
280  material.thermal_expansion_matrix(2);
281 
282  const libMesh::Point point(2.3, 3.1, 5.2);
283  const Real time = 2.34;
284  RealMatrixX dD_texp;
285  mat_texp.derivative(alpha, point, time, dD_texp);
286 
287  // Hard-coded in the true value of material thermal expansion matrix
288  RealMatrixX dD_texp_true = RealMatrixX::Zero(3,1);
289  dD_texp_true(0,0) = 1.0;
290  dD_texp_true(1,0) = 1.0;
291 
292  // Convert the test and truth Eigen::Matrix objects to std::vector
293  // since Catch2 has built in methods to compare vectors
294  std::vector<double> test = TEST::eigen_matrix_to_std_vector(dD_texp);
295  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(dD_texp_true);
296 
297  // Floating point approximations are diffcult to compare since the
298  // values typically aren't exactly equal due to numerical error.
299  // Therefore, we use the Approx comparison instead of Equals
300  CHECK_THAT( test, Catch::Approx<double>(truth) );
301  }
302 }
303 
304 
305 TEST_CASE("constant_isotropic_structural_material_2d_sensitivity",
306  "[isotropic],[material],[2D],[structural],[constant],[sensitivity]")
307 {
308  // Define Material Properties as MAST Parameters
309  MAST::Parameter E("E_param", 72.0e9); // Modulus of Elasticity
310  MAST::Parameter nu("nu_param", 0.33); // Poisson's ratio
311  MAST::Parameter kappa("kappa_param", 5.0/6.0); // Shear coefficient
312 
313  // Create field functions to dsitribute these constant parameters throughout the model
314  MAST::ConstantFieldFunction E_f("E", E);
315  MAST::ConstantFieldFunction nu_f("nu", nu);
316  MAST::ConstantFieldFunction kappa_f("kappa", kappa);
317 
318  // Initialize the material
320 
321  // Add the material property constant field functions to the material card
322  material.add(E_f);
323  material.add(nu_f);
324  material.add(kappa_f);
325 
326  // Ensure that the material properties were added before doing other checks
327  REQUIRE( material.contains("E") );
328  REQUIRE( material.contains("nu") );
329  REQUIRE( material.contains("kappa") );
330 
331  SECTION("material depends on the parameters that it should")
332  {
333  CHECK( material.depends_on(E) );
334  CHECK( material.depends_on(nu) );
335  CHECK( material.depends_on(kappa) );
336  }
337 
338  SECTION("material does not depend on other parameters")
339  {
340  MAST::Parameter dummy("dummy", 1.0);
341  CHECK_FALSE( material.depends_on(dummy) );
342  }
343 
344  SECTION("2D plane stress material constitutive matrix derivative w.r.t. elastic modulus is correct")
345  {
346  const bool is_plane_stress = true;
347  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
348  material.stiffness_matrix(2, is_plane_stress);
349 
350  const libMesh::Point point(2.3, 3.1, 5.2);
351  const Real time = 2.34;
352  RealMatrixX dK_mat;
353  mat_stiffness.derivative(E, point, time, dK_mat);
354 
355  // Hard-code in the true value of the material stiffness
356  RealMatrixX dK_mat_true = RealMatrixX::Zero(3,3);
357  dK_mat_true(0,0) = 1.122208506340478;
358  dK_mat_true(1,1) = 1.122208506340478;
359  dK_mat_true(0,1) = 0.370328807092358;
360  dK_mat_true(1,0) = 0.370328807092358;
361  dK_mat_true(2,2) = 0.375939849624060;
362 
363  // Convert the test and truth Eigen::Matrix objects to std::vector
364  // since Catch2 has built in methods to compare vectors
365  std::vector<double> test(dK_mat.data(), dK_mat.data()+dK_mat.rows()*dK_mat.cols());
366  std::vector<double> truth(dK_mat_true.data(), dK_mat_true.data()+dK_mat_true.rows()*dK_mat_true.cols());
367 
368  // Floating point approximations are diffcult to compare since the
369  // values typically aren't exactly equal due to numerical error.
370  // Therefore, we use the Approx comparison instead of Equals
371  CHECK_THAT( test, Catch::Approx<double>(truth) );
372  }
373 
374 
375  SECTION("2D plane stress material constitutive matrix derivative w.r.t. poisson's ratio is correct")
376  {
377  const bool is_plane_stress = true;
378  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
379  material.stiffness_matrix(2, is_plane_stress);
380 
381  const libMesh::Point point(2.3, 3.1, 5.2);
382  const Real time = 2.34;
383  RealMatrixX dK_mat;
384  mat_stiffness.derivative(nu, point, time, dK_mat);
385 
386  // Hard-code in the true value of the material stiffness
387  RealMatrixX dK_mat_true = RealMatrixX::Zero(3,3);
388  dK_mat_true(0,0) = 0.598444037945231e+11;
389  dK_mat_true(1,1) = 0.598444037945231e+11;
390  dK_mat_true(0,1) = 1.005476657087070e+11;
391  dK_mat_true(1,0) = 1.005476657087070e+11;
392  dK_mat_true(2,2) = -0.203516309570920e+11;
393 
394  // Convert the test and truth Eigen::Matrix objects to std::vector
395  // since Catch2 has built in methods to compare vectors
396  std::vector<double> test(dK_mat.data(), dK_mat.data()+dK_mat.rows()*dK_mat.cols());
397  std::vector<double> truth(dK_mat_true.data(), dK_mat_true.data()+dK_mat_true.rows()*dK_mat_true.cols());
398 
399  // Floating point approximations are diffcult to compare since the
400  // values typically aren't exactly equal due to numerical error.
401  // Therefore, we use the Approx comparison instead of Equals
402  CHECK_THAT( test, Catch::Approx<double>(truth) );
403  }
404 
405  SECTION("2D plane strain material constitutive matrix derivative w.r.t. elastic modulus is correct")
406  {
407  const bool is_plane_stress = false;
408  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
409  material.stiffness_matrix(2, is_plane_stress);
410 
411  const libMesh::Point point(2.3, 3.1, 5.2);
412  const Real time = 2.34;
413  RealMatrixX dK_mat;
414  mat_stiffness.derivative(E, point, time, dK_mat);
415 
416 // // 4th order central finite difference approximation
417 // RealMatrixX K_mat_h, K_mat_h2, K_mat_n, K_mat_n2;
418 // Real delta = 1e9;
419 // const Real E0 = E();
420 //
421 // E = E0 + delta;
422 // mat_stiffness.derivative(E, point, time, K_mat_h);
423 //
424 // E = E0 + 2.0*delta;
425 // mat_stiffness.derivative(E, point, time, K_mat_h2);
426 //
427 // E = E0 - delta;
428 // mat_stiffness.derivative(E, point, time, K_mat_n);
429 //
430 // E = E0 - 2.0*delta;
431 // mat_stiffness.derivative(E, point, time, K_mat_n2);
432 //
433 // E = E0;
434 //
435 // RealMatrixX dK_mat_true = (K_mat_n2 - 8.0*K_mat_n + 8.0*K_mat_h - K_mat_h2)/(12.0*delta);
436 
437  // Hard-code in the true value of the material stiffness
438  RealMatrixX dK_mat_true = RealMatrixX::Zero(3,3);
439  dK_mat_true(0,0) = 1.481645289694825;
440  dK_mat_true(1,1) = 1.481645289694825;
441  dK_mat_true(0,1) = 0.729765590446705;
442  dK_mat_true(1,0) = 0.729765590446705;
443  dK_mat_true(2,2) = 0.375939849624060;
444 
445  // Convert the test and truth Eigen::Matrix objects to std::vector
446  // since Catch2 has built in methods to compare vectors
447  std::vector<double> test(dK_mat.data(), dK_mat.data()+dK_mat.rows()*dK_mat.cols());
448  std::vector<double> truth(dK_mat_true.data(), dK_mat_true.data()+dK_mat_true.rows()*dK_mat_true.cols());
449 
450  // Floating point approximations are diffcult to compare since the
451  // values typically aren't exactly equal due to numerical error.
452  // Therefore, we use the Approx comparison instead of Equals
453  CHECK_THAT( test, Catch::Approx<double>(truth) );
454  }
455 
456  SECTION("2D plane strain material constitutive matrix derivative w.r.t. poisson's ratio is correct")
457  {
458  const bool is_plane_stress = false;
459  const MAST::FieldFunction<RealMatrixX>& mat_stiffness =
460  material.stiffness_matrix(2, is_plane_stress);
461 
462  const libMesh::Point point(2.3, 3.1, 5.2);
463  const Real time = 2.34;
464  RealMatrixX dK_mat;
465  mat_stiffness.derivative(nu, point, time, dK_mat);
466 
467  // Hard-code in the true value of the material stiffness
468  RealMatrixX dK_mat_true = RealMatrixX::Zero(3,3);
469  dK_mat_true(0,0) = 3.880894055520205e+11;
470  dK_mat_true(1,1) = 3.880894055520205e+11;
471  dK_mat_true(0,1) = 4.287926674662045e+11;
472  dK_mat_true(1,0) = 4.287926674662045e+11;
473  dK_mat_true(2,2) = -0.203516309570920e+11;
474 
475  // Convert the test and truth Eigen::Matrix objects to std::vector
476  // since Catch2 has built in methods to compare vectors
477  std::vector<double> test(dK_mat.data(), dK_mat.data()+dK_mat.rows()*dK_mat.cols());
478  std::vector<double> truth(dK_mat_true.data(), dK_mat_true.data()+dK_mat_true.rows()*dK_mat_true.cols());
479 
480  // Floating point approximations are diffcult to compare since the
481  // values typically aren't exactly equal due to numerical error.
482  // Therefore, we use the Approx comparison instead of Equals
483  CHECK_THAT( test, Catch::Approx<double>(truth) );
484  }
485 
486  SECTION("transverse shear material stiffness matrix derivative w.r.t. elastic modulus")
487  {
488  const MAST::FieldFunction<RealMatrixX>& mat_trans_shear =
490 
491  const libMesh::Point point(2.3, 3.1, 5.2);
492  const Real time = 2.34;
493  RealMatrixX dD_trans_shear;
494  mat_trans_shear.derivative(E, point, time, dD_trans_shear);
495 
496  // Hard-coded in the true value of the material transverse shear stiffness
497  RealMatrixX dD_trans_shear_true = RealMatrixX::Zero(2,2);
498  dD_trans_shear_true(0,0) = 0.313283208020050;
499  dD_trans_shear_true(1,1) = 0.313283208020050;
500 
501  // Convert the test and truth Eigen::Matrix objects to std::vector
502  // since Catch2 has built in methods to compare vectors
503  std::vector<double> test(dD_trans_shear.data(), dD_trans_shear.data()+dD_trans_shear.rows()*dD_trans_shear.cols());
504  std::vector<double> truth(dD_trans_shear_true.data(), dD_trans_shear_true.data()+dD_trans_shear_true.rows()*dD_trans_shear_true.cols());
505 
506  // Floating point approximations are diffcult to compare since the
507  // values typically aren't exactly equal due to numerical error.
508  // Therefore, we use the Approx comparison instead of Equals
509  CHECK_THAT( test, Catch::Approx<double>(truth) );
510  }
511 
512  SECTION("transverse shear material stiffness matrix derivative w.r.t. poisson's ratio")
513  {
514  const MAST::FieldFunction<RealMatrixX>& mat_trans_shear =
516 
517  const libMesh::Point point(2.3, 3.1, 5.2);
518  const Real time = 2.34;
519  RealMatrixX dD_trans_shear;
520  mat_trans_shear.derivative(nu, point, time, dD_trans_shear);
521 
522  // Hard-coded in the true value of the material transverse shear stiffness
523  RealMatrixX dD_trans_shear_true = RealMatrixX::Zero(2,2);
524  dD_trans_shear_true(0,0) = -1.695969246424331e+10;
525  dD_trans_shear_true(1,1) = -1.695969246424331e+10;
526 
527  // Convert the test and truth Eigen::Matrix objects to std::vector
528  // since Catch2 has built in methods to compare vectors
529  std::vector<double> test(dD_trans_shear.data(), dD_trans_shear.data()+dD_trans_shear.rows()*dD_trans_shear.cols());
530  std::vector<double> truth(dD_trans_shear_true.data(), dD_trans_shear_true.data()+dD_trans_shear_true.rows()*dD_trans_shear_true.cols());
531 
532  // Floating point approximations are diffcult to compare since the
533  // values typically aren't exactly equal due to numerical error.
534  // Therefore, we use the Approx comparison instead of Equals
535  CHECK_THAT( test, Catch::Approx<double>(truth) );
536  }
537 
538  SECTION("transverse shear material stiffness matrix derivative w.r.t. shear coefficient")
539  {
540  const MAST::FieldFunction<RealMatrixX>& mat_trans_shear =
542 
543  const libMesh::Point point(2.3, 3.1, 5.2);
544  const Real time = 2.34;
545  RealMatrixX dD_trans_shear;
546  mat_trans_shear.derivative(kappa, point, time, dD_trans_shear);
547 
548  // Hard-coded in the true value of the material transverse shear stiffness
549  RealMatrixX dD_trans_shear_true = RealMatrixX::Zero(2,2);
550  dD_trans_shear_true(0,0) = 2.706766917293233e+10;
551  dD_trans_shear_true(1,1) = 2.706766917293233e+10;
552 
553  // Convert the test and truth Eigen::Matrix objects to std::vector
554  // since Catch2 has built in methods to compare vectors
555  std::vector<double> test = TEST::eigen_matrix_to_std_vector(dD_trans_shear);
556  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(dD_trans_shear_true);
557 
558  // Floating point approximations are diffcult to compare since the
559  // values typically aren't exactly equal due to numerical error.
560  // Therefore, we use the Approx comparison instead of Equals
561  CHECK_THAT( test, Catch::Approx<double>(truth) );
562  }
563 }
564 
565 
566 TEST_CASE("constant_isotropic_transient_heat_transfer_material_2d_sensitivity",
567  "[isotropic],[material],[2D],[heat_transfer],[constant],[sensitivity],[transient]")
568 {
569  // Define Material Properties as MAST Parameters
570  MAST::Parameter rho("rho_param", 1234.5); // Density
571  MAST::Parameter cp("cp_param", 908.0); // Specific Heat Capacity
572  MAST::Parameter k("k_param", 237.0); // Thermal Conductivity
573 
574  // Create field functions to dsitribute these constant parameters throughout the model
575  MAST::ConstantFieldFunction rho_f("rho", rho);
576  MAST::ConstantFieldFunction cp_f("cp", cp);
577  MAST::ConstantFieldFunction k_f("k_th", k);
578 
579  // Initialize the material
581 
582  // Add the material property constant field functions to the material card
583  material.add(rho_f);
584  material.add(cp_f);
585  material.add(k_f);
586 
587  // Ensure that the material properties were added before doing other checks
588  REQUIRE( material.contains("rho") );
589  REQUIRE( material.contains("cp") );
590  REQUIRE( material.contains("k_th") );
591 
592  SECTION("material depends on the parameters that it should")
593  {
594  CHECK( material.depends_on(rho) );
595  CHECK( material.depends_on(cp) );
596  CHECK( material.depends_on(k) );
597  }
598 
599  SECTION("material does not depend on other parameters")
600  {
601  MAST::Parameter dummy("dummy", 1.0);
602  CHECK_FALSE( material.depends_on(dummy) );
603  }
604 
605  SECTION("2D thermal capacitance matrix derivative w.r.t. specifc heat capacity")
606  {
607  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
608  material.capacitance_matrix(2);
609 
610  const libMesh::Point point(2.3, 3.1, 5.2);
611  const Real time = 2.34;
612  RealMatrixX dD_cp;
613  mat_capacit.derivative(cp, point, time, dD_cp);
614 
615  // Hard-coded values for thermal conductivity matrix
616  RealMatrixX dD_cp_true = RealMatrixX::Identity(1,1);
617  dD_cp_true *= 1234.5;
618 
619  // Convert the test and truth Eigen::Matrix objects to std::vector
620  // since Catch2 has built in methods to compare vectors
621  std::vector<double> test = TEST::eigen_matrix_to_std_vector(dD_cp);
622  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(dD_cp_true);
623 
624  // Floating point approximations are diffcult to compare since the
625  // values typically aren't exactly equal due to numerical error.
626  // Therefore, we use the Approx comparison instead of Equals
627  CHECK_THAT( test, Catch::Approx<double>(truth) );
628  }
629 
630  SECTION("2D thermal capacitance matrix derivative w.r.t. density")
631  {
632  const MAST::FieldFunction<RealMatrixX>& mat_capacit =
633  material.capacitance_matrix(2);
634 
635  const libMesh::Point point(2.3, 3.1, 5.2);
636  const Real time = 2.34;
637  RealMatrixX dD_cp;
638  mat_capacit.derivative(rho, point, time, dD_cp);
639 
640  // Hard-coded values for thermal conductivity matrix
641  RealMatrixX dD_cp_true = RealMatrixX::Identity(1,1);
642  dD_cp_true *= 908.0;
643 
644  // Convert the test and truth Eigen::Matrix objects to std::vector
645  // since Catch2 has built in methods to compare vectors
646  std::vector<double> test = TEST::eigen_matrix_to_std_vector(dD_cp);
647  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(dD_cp_true);
648 
649  // Floating point approximations are diffcult to compare since the
650  // values typically aren't exactly equal due to numerical error.
651  // Therefore, we use the Approx comparison instead of Equals
652  CHECK_THAT( test, Catch::Approx<double>(truth) );
653  }
654 }
655 
656 
657 TEST_CASE("constant_isotropic_heat_transfer_material_2d_sensitivity",
658  "[isotropic],[material],[2D],[heat_transfer][constant],[sensitivity]")
659 {
660  // Define Material Properties as MAST Parameters
661  MAST::Parameter rho("rho_param", 1234.5); // Density
662  MAST::Parameter cp("cp_param", 908.0); // Specific Heat Capacity
663  MAST::Parameter k("k_param", 237.0); // Thermal Conductivity
664 
665  // Create field functions to dsitribute these constant parameters throughout the model
666  MAST::ConstantFieldFunction rho_f("rho", rho);
667  MAST::ConstantFieldFunction cp_f("cp", cp);
668  MAST::ConstantFieldFunction k_f("k_th", k);
669 
670  // Initialize the material
672 
673  // Add the material property constant field functions to the material card
674  material.add(rho_f);
675  material.add(cp_f);
676  material.add(k_f);
677 
678  // Ensure that the material properties were added before doing other checks
679  REQUIRE( material.contains("rho") );
680  REQUIRE( material.contains("cp") );
681  REQUIRE( material.contains("k_th") );
682 
683  SECTION("material depends on the parameters that it should")
684  {
685  CHECK( material.depends_on(rho) );
686  CHECK( material.depends_on(cp) );
687  CHECK( material.depends_on(k) );
688  }
689 
690  SECTION("material does not depend on other parameters")
691  {
692  MAST::Parameter dummy("dummy", 1.0);
693  CHECK_FALSE( material.depends_on(dummy) );
694  }
695 
696  SECTION("2D thermal conductivity matrix derivative w.r.t. thermal conductivity")
697  {
698  const MAST::FieldFunction<RealMatrixX>& mat_conduct =
699  material.conductance_matrix(2);
700 
701  const libMesh::Point point(2.3, 3.1, 5.2);
702  const Real time = 2.34;
703  RealMatrixX dD_k;
704  mat_conduct.derivative(k, point, time, dD_k);
705 
706  // Hard-coded values for thermal conductivity matrix
707  RealMatrixX dD_k_true = RealMatrixX::Identity(2,2);
708 
709  // Convert the test and truth Eigen::Matrix objects to std::vector
710  // since Catch2 has built in methods to compare vectors
711  std::vector<double> test = TEST::eigen_matrix_to_std_vector(dD_k);
712  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(dD_k_true);
713 
714  // Floating point approximations are diffcult to compare since the
715  // values typically aren't exactly equal due to numerical error.
716  // Therefore, we use the Approx comparison instead of Equals
717  CHECK_THAT( test, Catch::Approx<double>(truth) );
718  }
719 }
TEST_CASE("constant_isotropic_thermoelastic_material_1d_sensitivity", "[isotropic],[material],[constant],[1D],[thermoelastic],[sensitivity]")
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 void derivative(const MAST::FunctionBase &f, ValType &v) const
calculates the value of the function derivative and returns it in v.
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)
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)
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