MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_solid_2d_section_element_property_card.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"
12 
13 // Custom includes
14 #include "test_helpers.h"
15 
16 // TODO: Need to test with other types of materials.
17 // TODO: Need to test function that gets material.
18 
19 extern libMesh::LibMeshInit* p_global_init;
20 
21 TEST_CASE("element_property_card_constant_heat_transfer_2d",
22  "[heat_transfer],[2D],[isotropic],[constant],[property]")
23 {
24  const uint dim = 2;
25 
26  // Define Material Properties as MAST Parameters
27  MAST::Parameter k("k_param", 237.0); // Thermal Conductivity
28 
29  // Define Section Properties as MAST Parameters
30  MAST::Parameter thickness("th_param", 0.06); // Section thickness
31  MAST::Parameter offset("off_param", 0.03); // Section offset
32 
33  // Create field functions to dsitribute these constant parameters throughout the model
34  MAST::ConstantFieldFunction k_f("k_th", k);
35  MAST::ConstantFieldFunction thickness_f("h", thickness);
36  MAST::ConstantFieldFunction offset_f("off", offset);
37 
38  // Initialize the material
40 
41  // Add the material property constant field functions to the material card
42  material.add(k_f);
43 
44  // Initialize the section
46 
47  // Add the section property constant field functions to the section card
48  section.add(thickness_f);
49  section.add(offset_f);
50 
51  // Add the material card to the section card
52  section.set_material(material);
53 
54  REQUIRE( section.dim() == dim); // Ensure section is 2 dimensional
55  REQUIRE( section.depends_on(thickness) );
56  REQUIRE( section.depends_on(offset) );
57  REQUIRE( section.depends_on(k) );
58 
59  SECTION("2D section thermal conductance matrix")
60  {
71  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> conduct_mat = section.thermal_conductance_matrix();
72 
73  const libMesh::Point point(2.3, 3.1, 5.2);
74  const Real time = 2.34;
75  RealMatrixX D_sec_conduc;
76  conduct_mat->operator()(point, time, D_sec_conduc);
77 
78  // Hard-coded value of the section's extension stiffness
79  RealMatrixX D_sec_conduc_true = RealMatrixX::Zero(2,2);
80  D_sec_conduc_true(0,0) = 237.0*0.06;
81  D_sec_conduc_true(1,1) = 237.0*0.06;
82 
83  // Convert the test and truth Eigen::Matrix objects to std::vector
84  // since Catch2 has built in methods to compare vectors
85  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_conduc);
86  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_conduc_true);
87 
88  // Floating point approximations are diffcult to compare since the
89  // values typically aren't exactly equal due to numerical error.
90  // Therefore, we use the Approx comparison instead of Equals
91  CHECK_THAT( test, Catch::Approx<double>(truth) );
92  }
93 }
94 
95 
96 TEST_CASE("element_property_card_constant_transient_heat_transfer_2d",
97  "[heat_transfer],[2D],[isotropic],[constant],[property],[transient]")
98 {
99  const uint dim = 2;
100 
101  // Define Material Properties as MAST Parameters
102  MAST::Parameter rho("rho_param", 1420.5); // Density
103  MAST::Parameter cp("cp_param", 908.0); // Specific Heat Capacity
104 
105  // Define Section Properties as MAST Parameters
106  MAST::Parameter thickness("th_param", 0.06); // Section thickness
107  MAST::Parameter offset("off_param", 0.03); // Section offset
108 
109  // Create field functions to dsitribute these constant parameters throughout the model
110  MAST::ConstantFieldFunction rho_f("rho", rho);
111  MAST::ConstantFieldFunction cp_f("cp", cp);
112  MAST::ConstantFieldFunction thickness_f("h", thickness);
113  MAST::ConstantFieldFunction offset_f("off", offset);
114 
115  // Initialize the material
117 
118  // Add the material property constant field functions to the material card
119  material.add(rho_f);
120  material.add(cp_f);
121 
122  // Initialize the section
124 
125  // Add the section property constant field functions to the section card
126  section.add(thickness_f);
127  section.add(offset_f);
128 
129  // Add the material card to the section card
130  section.set_material(material);
131 
132  REQUIRE( section.dim() == dim); // Ensure section is 2 dimensional
133  REQUIRE( section.depends_on(thickness) );
134  REQUIRE( section.depends_on(offset) );
135  REQUIRE( section.depends_on(cp) );
136  REQUIRE( section.depends_on(rho) );
137 
138  SECTION("2D section thermal capacitance matrix")
139  {
150  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> capaci_mat = section.thermal_capacitance_matrix();
151 
152  const libMesh::Point point(2.3, 3.1, 5.2);
153  const Real time = 2.34;
154  RealMatrixX D_sec_capac;
155  capaci_mat->operator()(point, time, D_sec_capac);
156 
157  // Hard-coded value of the section's extension stiffness
158  RealMatrixX D_sec_capac_true = RealMatrixX::Zero(1,1);
159  D_sec_capac_true(0,0) = 908.0*1420.5*0.06;
160 
161  // Convert the test and truth Eigen::Matrix objects to std::vector
162  // since Catch2 has built in methods to compare vectors
163  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_capac);
164  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_capac_true);
165 
166  // Floating point approximations are diffcult to compare since the
167  // values typically aren't exactly equal due to numerical error.
168  // Therefore, we use the Approx comparison instead of Equals
169  CHECK_THAT( test, Catch::Approx<double>(truth) );
170  }
171 }
172 
173 
174 TEST_CASE("element_property_card_constant_thermoelastic_2d",
175  "[thermoelastic],[2D],[isotropic],[constant],[property]")
176 {
177  const uint dim = 2;
178 
179  // Define Material Properties as MAST Parameters
180  MAST::Parameter E("E_param", 72.0e9); // Modulus of Elasticity
181  MAST::Parameter nu("nu_param", 0.33); // Poisson's ratio
182  MAST::Parameter alpha("alpha_param", 5.43e-05); // Coefficient of thermal expansion
183 
184  // Define Section Properties as MAST Parameters
185  MAST::Parameter thickness("th_param", 0.06); // Section thickness
186  MAST::Parameter offset("off_param", 0.03); // Section offset
187 
188  // Create field functions to dsitribute these constant parameters throughout the model
189  MAST::ConstantFieldFunction E_f("E", E);
190  MAST::ConstantFieldFunction nu_f("nu", nu);
191  MAST::ConstantFieldFunction alpha_f("alpha_expansion", alpha);
192  MAST::ConstantFieldFunction thickness_f("h", thickness);
193  MAST::ConstantFieldFunction offset_f("off", offset);
194 
195  // Initialize the material
197 
198  // Add the material property constant field functions to the material card
199  material.add(E_f);
200  material.add(nu_f);
201  material.add(alpha_f);
202 
203  // Initialize the section
205 
206  // Add the section property constant field functions to the section card
207  section.add(thickness_f);
208  section.add(offset_f);
209 
210  // Add the material card to the section card
211  section.set_material(material);
212 
213  REQUIRE( section.dim() == dim); // Ensure section is 2 dimensional
214  REQUIRE( section.depends_on(thickness) );
215  REQUIRE( section.depends_on(offset) );
216  REQUIRE( section.depends_on(alpha) );
217 
218  SECTION("2D plane stress thermal expansion A matrix")
219  {
224  section.set_plane_stress(true);
225  REQUIRE( section.plane_stress() );
226 
237  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> texp_A_mat = section.thermal_expansion_A_matrix();
238 
239  const libMesh::Point point(2.3, 3.1, 5.2);
240  const Real time = 2.34;
241  RealMatrixX D_sec_texpA;
242  texp_A_mat->operator()(point, time, D_sec_texpA);
243 
244  // Hard-coded value of the section's extension stiffness
245  RealMatrixX D_sec_texpA_true = RealMatrixX::Zero(3,1);
246  D_sec_texpA_true(0,0) = 3.501134328358209e+05;
247  D_sec_texpA_true(1,0) = 3.501134328358209e+05;
248 
249  // Convert the test and truth Eigen::Matrix objects to std::vector
250  // since Catch2 has built in methods to compare vectors
251  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_texpA);
252  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_texpA_true);
253 
254  // Floating point approximations are diffcult to compare since the
255  // values typically aren't exactly equal due to numerical error.
256  // Therefore, we use the Approx comparison instead of Equals
257  CHECK_THAT( test, Catch::Approx<double>(truth) );
258  }
259 
260 
261  SECTION("2D plane strain thermal expansion A matrix")
262  {
267  section.set_plane_stress(false);
268  REQUIRE_FALSE( section.plane_stress() );
269 
280  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> texp_A_mat = section.thermal_expansion_A_matrix();
281 
282  const libMesh::Point point(2.3, 3.1, 5.2);
283  const Real time = 2.34;
284  RealMatrixX D_sec_texpA;
285  texp_A_mat->operator()(point, time, D_sec_texpA);
286 
287  // Hard-coded value of the section's extension stiffness
288  RealMatrixX D_sec_texpA_true = RealMatrixX::Zero(3,1);
289  D_sec_texpA_true(0,0) = 5.187439186200796e+05;
290  D_sec_texpA_true(1,0) = 5.187439186200796e+05;
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(D_sec_texpA);
295  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_texpA_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  SECTION("2D plane stress thermal expansion B matrix")
305  {
310  section.set_plane_stress(true);
311  REQUIRE( section.plane_stress() );
312 
323  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> texp_B_mat = section.thermal_expansion_B_matrix();
324 
325  const libMesh::Point point(2.3, 3.1, 5.2);
326  const Real time = 2.34;
327  RealMatrixX D_sec_texpB;
328  texp_B_mat->operator()(point, time, D_sec_texpB);
329 
330  // Hard-coded value of the section's extension stiffness
331  RealMatrixX D_sec_texpB_true = RealMatrixX::Zero(3,1);
332  D_sec_texpB_true(0,0) = 1.050340298507463e+04;
333  D_sec_texpB_true(1,0) = 1.050340298507463e+04;
334 
335  // Convert the test and truth Eigen::Matrix objects to std::vector
336  // since Catch2 has built in methods to compare vectors
337  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_texpB);
338  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_texpB_true);
339 
340  // Floating point approximations are diffcult to compare since the
341  // values typically aren't exactly equal due to numerical error.
342  // Therefore, we use the Approx comparison instead of Equals
343  CHECK_THAT( test, Catch::Approx<double>(truth) );
344  }
345 
346  SECTION("2D plane strain thermal expansion B matrix")
347  {
352  section.set_plane_stress(false);
353  REQUIRE_FALSE( section.plane_stress() );
354 
365  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> texp_B_mat = section.thermal_expansion_B_matrix();
366 
367  const libMesh::Point point(2.3, 3.1, 5.2);
368  const Real time = 2.34;
369  RealMatrixX D_sec_texpB;
370  texp_B_mat->operator()(point, time, D_sec_texpB);
371 
372  // Hard-coded value of the section's extension stiffness
373  RealMatrixX D_sec_texpB_true = RealMatrixX::Zero(3,1);
374  D_sec_texpB_true(0,0) = 1.556231755860239e+04;
375  D_sec_texpB_true(1,0) = 1.556231755860239e+04;
376 
377  // Convert the test and truth Eigen::Matrix objects to std::vector
378  // since Catch2 has built in methods to compare vectors
379  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_texpB);
380  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_texpB_true);
381 
382  // Floating point approximations are diffcult to compare since the
383  // values typically aren't exactly equal due to numerical error.
384  // Therefore, we use the Approx comparison instead of Equals
385  CHECK_THAT( test, Catch::Approx<double>(truth) );
386  }
387 }
388 
389 
390 TEST_CASE("element_property_card_constant_dynamic_2d",
391  "[dynamic],[2D],[isotropic],[constant],[property]")
392 {
393  const uint dim = 2;
394 
395  // Define Material Properties as MAST Parameters
396  MAST::Parameter rho("rho_param", 1420.5); // Density
397  MAST::Parameter E("E_param", 72.0e9); // Modulus of Elasticity
398  MAST::Parameter nu("nu_param", 0.33); // Poisson's ratio
399  MAST::Parameter kappa("kappa_param", 5.0/6.0); // Shear coefficient
400 
401  // Define Section Properties as MAST Parameters
402  MAST::Parameter thickness("th_param", 0.06); // Section thickness
403  MAST::Parameter offset("off_param", 0.03); // Section offset
404 
405  // Create field functions to dsitribute these constant parameters throughout the model
406  MAST::ConstantFieldFunction rho_f("rho", rho);
407  MAST::ConstantFieldFunction E_f("E", E);
408  MAST::ConstantFieldFunction nu_f("nu", nu);
409  MAST::ConstantFieldFunction kappa_f("kappa", kappa);
410  MAST::ConstantFieldFunction thickness_f("h", thickness);
411  MAST::ConstantFieldFunction offset_f("off", offset);
412 
413  // Initialize the material
415 
416  // Add the material property constant field functions to the material card
417  material.add(rho_f);
418  material.add(E_f);
419  material.add(nu_f);
420 
421  // Initialize the section
423 
424  // Add the section property constant field functions to the section card
425  section.add(thickness_f);
426  section.add(offset_f);
427  section.add(kappa_f);
428 
429  // Add the material card to the section card
430  section.set_material(material);
431 
432  REQUIRE( section.dim() == dim); // Ensure section is 2 dimensional
433  REQUIRE( section.depends_on(thickness) );
434  REQUIRE( section.depends_on(offset) );
435  REQUIRE( section.depends_on(E) );
436  REQUIRE( section.depends_on(nu) );
437  REQUIRE( section.depends_on(kappa) );
438 
439  REQUIRE_FALSE( section.if_diagonal_mass_matrix() );
440 
441  section.set_diagonal_mass_matrix(true);
442  REQUIRE( section.if_diagonal_mass_matrix() );
443 
444  SECTION("2D section inertia matrix")
445  {
456  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> inertia_mat = section.inertia_matrix();
457 
458  const libMesh::Point point(2.3, 3.1, 5.2);
459  const Real time = 2.34;
460  RealMatrixX D_sec_iner;
461  inertia_mat->operator()(point, time, D_sec_iner);
462 
463  // Hard-coded value of the section's extension stiffness
464  RealMatrixX D_sec_iner_true = RealMatrixX::Zero(6,6);
465  D_sec_iner_true(0,0) = D_sec_iner_true(1,1) = D_sec_iner_true(2,2) = 0.06;
466  D_sec_iner_true(0,4) = D_sec_iner_true(4,0) = 0.03*0.06;
467  D_sec_iner_true(1,3) = D_sec_iner_true(3,1) = -0.03*0.06;
468  D_sec_iner_true(3,3) = D_sec_iner_true(4,4) = pow(0.06,3.0)/12.0 + 0.06*0.03*0.03;
469  D_sec_iner_true(5,5) = pow(0.06,3.0)/12.0 * 1.0e-6; // Ignore this term
470 
471  D_sec_iner_true *= 1420.5;
472 
473  // Convert the test and truth Eigen::Matrix objects to std::vector
474  // since Catch2 has built in methods to compare vectors
475  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_iner);
476  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_iner_true);
477 
478  // Floating point approximations are diffcult to compare since the
479  // values typically aren't exactly equal due to numerical error.
480  // Therefore, we use the Approx comparison instead of Equals
481  CHECK_THAT( test, Catch::Approx<double>(truth) );
482  }
483 }
484 
485 
486 
487 TEST_CASE("element_property_card_constant_structural_2d",
488  "[structural],[2D],[isotropic],[constant],[property]")
489 {
490  const uint dim = 2;
491 
492  // Define Material Properties as MAST Parameters
493  MAST::Parameter E("E_param", 72.0e9); // Modulus of Elasticity
494  MAST::Parameter nu("nu_param", 0.33); // Poisson's ratio
495  MAST::Parameter kappa("kappa_param", 5.0/6.0); // Shear coefficient
496 
497  // Define Section Properties as MAST Parameters
498  MAST::Parameter thickness("th_param", 0.06); // Section thickness
499  MAST::Parameter offset("off_param", 0.03); // Section offset
500 
501  // Create field functions to dsitribute these constant parameters throughout the model
502  MAST::ConstantFieldFunction E_f("E", E);
503  MAST::ConstantFieldFunction nu_f("nu", nu);
504  MAST::ConstantFieldFunction kappa_f("kappa", kappa);
505  MAST::ConstantFieldFunction thickness_f("h", thickness);
506  MAST::ConstantFieldFunction offset_f("off", offset);
507 
508  // Initialize the material
510 
511  // Add the material property constant field functions to the material card
512  material.add(E_f);
513  material.add(nu_f);
514 
515  // Initialize the section
517 
518  // Add the section property constant field functions to the section card
519  section.add(thickness_f);
520  section.add(offset_f);
521  section.add(kappa_f);
522 
523  // Add the material card to the section card
524  section.set_material(material);
525 
526  REQUIRE( section.dim() == dim); // Ensure section is 2 dimensional
527  REQUIRE( section.depends_on(thickness) );
528  REQUIRE( section.depends_on(offset) );
529  REQUIRE( section.depends_on(E) );
530  REQUIRE( section.depends_on(nu) );
531  REQUIRE( section.depends_on(kappa) );
532 
533  SECTION("solid_2d_section is isotropic")
534  {
535  CHECK( section.if_isotropic() );
536  }
537 
538  SECTION("set_get_strain_type")
539  {
540  // Check that the default is linear strain
541  REQUIRE( section.strain_type() == MAST::LINEAR_STRAIN );
542 
544  REQUIRE( section.strain_type() == MAST::LINEAR_STRAIN );
545  REQUIRE_FALSE( section.strain_type() == MAST::NONLINEAR_STRAIN );
546 
548  REQUIRE( section.strain_type() == MAST::NONLINEAR_STRAIN );
549  REQUIRE_FALSE( section.strain_type() == MAST::LINEAR_STRAIN );
550  }
551 
552  SECTION("set_get_bending_model")
553  {
554  // NOTE: MAST::BERNOULLI AND MAST::TIMOSHENKO are not valid options for 2D sections, even though their input is accepted.
555  section.set_bending_model(MAST::DKT);
559 
560  // TODO: Implement element creation for testing of getting bending_model and checking default
561  //REQUIRE( section.bending_model()
562  }
563 
564  SECTION("set_get_plane_stress_strain")
565  {
566  // Check that default is plane stress
567  REQUIRE( section.plane_stress() == true );
568 
569  section.set_plane_stress(false);
570  REQUIRE( section.plane_stress() == false );
571 
572  section.set_plane_stress(true);
573  REQUIRE( section.plane_stress() == true);
574  }
575 
576 
577 // SECTION("quadrature_order")
578 // {
579 // section.set_bending_model(MAST::MINDLIN);
580 // REQUIRE( section.extra_quadrature_order(elem) == 0 );
581 //
582 // section.set_bending_model(MAST::DKT);
583 // REQUIRE( section.extra_quadrature_order(elem) == 2 );
584 // }
585 
586  SECTION("2D plane stress extension stiffness matrix")
587  {
592  section.set_plane_stress(true);
593  REQUIRE( section.plane_stress() );
594 
605  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> extension_stiffness_mat = section.stiffness_A_matrix();
606 
607  const libMesh::Point point(2.3, 3.1, 5.2);
608  const Real time = 2.34;
609  RealMatrixX D_sec_ext;
610  extension_stiffness_mat->operator()(point, time, D_sec_ext);
611 
612  // Hard-coded value of the section's extension stiffness
613  RealMatrixX D_sec_ext_true = RealMatrixX::Zero(3,3);
614  D_sec_ext_true(0,0) = 4.847940747390865e+09;
615  D_sec_ext_true(1,1) = 4.847940747390865e+09;
616  D_sec_ext_true(0,1) = 1.599820446638986e+09;
617  D_sec_ext_true(1,0) = 1.599820446638986e+09;
618  D_sec_ext_true(2,2) = 1.624060150375940e+09;
619 
620  // Convert the test and truth Eigen::Matrix objects to std::vector
621  // since Catch2 has built in methods to compare vectors
622  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_ext);
623  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_ext_true);
624 
625  // Floating point approximations are diffcult to compare since the
626  // values typically aren't exactly equal due to numerical error.
627  // Therefore, we use the Approx comparison instead of Equals
628  CHECK_THAT( test, Catch::Approx<double>(truth) );
629  }
630 
631 
632  SECTION("2D plane strain extension stiffness matrix")
633  {
638  section.set_plane_stress(false);
639  REQUIRE_FALSE( section.plane_stress() );
640 
651  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> extension_stiffness_mat = section.stiffness_A_matrix();
652 
653  const libMesh::Point point(2.3, 3.1, 5.2);
654  const Real time = 2.34;
655  RealMatrixX D_sec_ext;
656  extension_stiffness_mat->operator()(point, time, D_sec_ext);
657 
658  // Hard-coded value of the section's extension stiffness
659  RealMatrixX D_sec_ext_true = RealMatrixX::Zero(3,3);
660  D_sec_ext_true(0,0) = 6.400707651481644e+09;
661  D_sec_ext_true(1,1) = 6.400707651481644e+09;
662  D_sec_ext_true(0,1) = 3.152587350729766e+09;
663  D_sec_ext_true(1,0) = 3.152587350729766e+09;
664  D_sec_ext_true(2,2) = 1.624060150375940e+09;
665 
666  // Convert the test and truth Eigen::Matrix objects to std::vector
667  // since Catch2 has built in methods to compare vectors
668  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_ext);
669  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_ext_true);
670 
671  // Floating point approximations are diffcult to compare since the
672  // values typically aren't exactly equal due to numerical error.
673  // Therefore, we use the Approx comparison instead of Equals
674  CHECK_THAT( test, Catch::Approx<double>(truth) );
675  }
676 
677 
678  SECTION("2D plane stress bending section stiffness matrix")
679  {
684  section.set_plane_stress(true);
685  REQUIRE( section.plane_stress() );
686 
697  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> bending_stiffness_mat = section.stiffness_D_matrix();
698 
699  const libMesh::Point point(2.3, 3.1, 5.2);
700  const Real time = 2.34;
701  RealMatrixX D_sec_bnd;
702  bending_stiffness_mat->operator()(point, time, D_sec_bnd);
703 
704  // Hard-coded value of the section's extension stiffness
705  RealMatrixX D_sec_bnd_true = RealMatrixX::Zero(3,3);
706  D_sec_bnd_true(0,0) = 5.817528896869037e+06;
707  D_sec_bnd_true(1,1) = 5.817528896869037e+06;
708  D_sec_bnd_true(0,1) = 1.919784535966782e+06;
709  D_sec_bnd_true(1,0) = 1.919784535966782e+06;
710  D_sec_bnd_true(2,2) = 1.948872180451127e+06;
711 
712  // Convert the test and truth Eigen::Matrix objects to std::vector
713  // since Catch2 has built in methods to compare vectors
714  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_bnd);
715  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_bnd_true);
716 
717  // Floating point approximations are diffcult to compare since the
718  // values typically aren't exactly equal due to numerical error.
719  // Therefore, we use the Approx comparison instead of Equals
720  CHECK_THAT( test, Catch::Approx<double>(truth) );
721  }
722 
723 
724  SECTION("2D plane strain bending section stiffness matrix")
725  {
730  section.set_plane_stress(false);
731  REQUIRE_FALSE( section.plane_stress() );
732 
743  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> bending_stiffness_mat = section.stiffness_D_matrix();
744 
745  const libMesh::Point point(2.3, 3.1, 5.2);
746  const Real time = 2.34;
747  RealMatrixX D_sec_bnd;
748  bending_stiffness_mat->operator()(point, time, D_sec_bnd);
749 
750  // Hard-coded value of the section's extension stiffness
751  RealMatrixX D_sec_bnd_true = RealMatrixX::Zero(3,3);
752  D_sec_bnd_true(0,0) = 7.680849181777972e+06;
753  D_sec_bnd_true(1,1) = 7.680849181777972e+06;
754  D_sec_bnd_true(0,1) = 3.783104820875719e+06;
755  D_sec_bnd_true(1,0) = 3.783104820875719e+06;
756  D_sec_bnd_true(2,2) = 1.948872180451128e+06;
757 
758  // Convert the test and truth Eigen::Matrix objects to std::vector
759  // since Catch2 has built in methods to compare vectors
760  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_bnd);
761  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_bnd_true);
762 
763  // Floating point approximations are diffcult to compare since the
764  // values typically aren't exactly equal due to numerical error.
765  // Therefore, we use the Approx comparison instead of Equals
766  CHECK_THAT( test, Catch::Approx<double>(truth) );
767  }
768 
769 
770  SECTION("2D plane stress extension-bending section stiffness matrix")
771  {
777  section.set_plane_stress(true);
778  REQUIRE( section.plane_stress() );
779 
790  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> bndext_stiffness_mat = section.stiffness_B_matrix();
791 
792  const libMesh::Point point(2.3, 3.1, 5.2);
793  const Real time = 2.34;
794  RealMatrixX D_sec_bndext;
795  bndext_stiffness_mat->operator()(point, time, D_sec_bndext);
796 
797  // Hard-coded value of the section's extension stiffness
798  RealMatrixX D_sec_bndext_true = RealMatrixX::Zero(3,3);
799  D_sec_bndext_true(0,0) = 1.454382224217259e+08;
800  D_sec_bndext_true(1,1) = 1.454382224217259e+08;
801  D_sec_bndext_true(0,1) = 0.479946133991696e+08;
802  D_sec_bndext_true(1,0) = 0.479946133991696e+08;
803  D_sec_bndext_true(2,2) = 0.487218045112782e+08;
804 
805  // Convert the test and truth Eigen::Matrix objects to std::vector
806  // since Catch2 has built in methods to compare vectors
807  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_bndext);
808  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_bndext_true);
809 
810  // Floating point approximations are diffcult to compare since the
811  // values typically aren't exactly equal due to numerical error.
812  // Therefore, we use the Approx comparison instead of Equals
813  CHECK_THAT( test, Catch::Approx<double>(truth) );
814  }
815 
816 
817  SECTION("2D plane strain extension-bending section stiffness matrix")
818  {
824  section.set_plane_stress(false);
825  REQUIRE_FALSE( section.plane_stress() );
826 
837  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> bndext_stiffness_mat = section.stiffness_B_matrix();
838 
839  const libMesh::Point point(2.3, 3.1, 5.2);
840  const Real time = 2.34;
841  RealMatrixX D_sec_bndext;
842  bndext_stiffness_mat->operator()(point, time, D_sec_bndext);
843 
844  // Hard-coded value of the section's extension stiffness
845  RealMatrixX D_sec_bndext_true = RealMatrixX::Zero(3,3);
846  D_sec_bndext_true(0,0) = 1.920212295444493e+08;
847  D_sec_bndext_true(1,1) = 1.920212295444493e+08;
848  D_sec_bndext_true(0,1) = 0.945776205218930e+08;
849  D_sec_bndext_true(1,0) = 0.945776205218930e+08;
850  D_sec_bndext_true(2,2) = 0.487218045112782e+08;
851 
852  // Convert the test and truth Eigen::Matrix objects to std::vector
853  // since Catch2 has built in methods to compare vectors
854  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_bndext);
855  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_bndext_true);
856 
857  // Floating point approximations are diffcult to compare since the
858  // values typically aren't exactly equal due to numerical error.
859  // Therefore, we use the Approx comparison instead of Equals
860  CHECK_THAT( test, Catch::Approx<double>(truth) );
861  }
862 
863 
864  SECTION("2D transverse shear section stiffness matrix")
865  {
871  section.set_plane_stress(true);
872  REQUIRE( section.plane_stress() );
873 
884  std::unique_ptr<MAST::FieldFunction<RealMatrixX>> trans_shear_stiffness_mat = section.transverse_shear_stiffness_matrix();
885 
886  const libMesh::Point point(2.3, 3.1, 5.2);
887  const Real time = 2.34;
888  RealMatrixX D_sec_shr;
889  trans_shear_stiffness_mat->operator()(point, time, D_sec_shr);
890 
891  // Hard-coded value of the section's extension stiffness
892  RealMatrixX D_sec_shr_true = RealMatrixX::Zero(2,2);
893  D_sec_shr_true(0,0) = 1.353383458646616e+09;
894  D_sec_shr_true(1,1) = 1.353383458646616e+09;
895 
896  // Convert the test and truth Eigen::Matrix objects to std::vector
897  // since Catch2 has built in methods to compare vectors
898  std::vector<double> test = TEST::eigen_matrix_to_std_vector(D_sec_shr);
899  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(D_sec_shr_true);
900 
901  // Floating point approximations are diffcult to compare since the
902  // values typically aren't exactly equal due to numerical error.
903  // Therefore, we use the Approx comparison instead of Equals
904  CHECK_THAT( test, Catch::Approx<double>(truth) );
905  }
906 }
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > stiffness_B_matrix(const MAST::ElementBase &e) const
void set_diagonal_mass_matrix(bool m)
sets the mass matrix to be diagonal or consistent
TEST_CASE("element_property_card_constant_heat_transfer_2d", "[heat_transfer],[2D],[isotropic],[constant],[property]")
virtual bool if_isotropic() const
return true if the property is isotropic
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > transverse_shear_stiffness_matrix(const MAST::ElementBase &e) const
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
void set_plane_stress(bool val)
sets the flag for plane stress.
libMesh::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26
libMesh::Real Real
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > stiffness_D_matrix(const MAST::ElementBase &e) const
virtual unsigned int dim() const
dimension of the element for which this property is defined
void set_strain(MAST::StrainType strain)
sets the type of strain to be used, which is LINEAR_STRAIN by default
void add(MAST::FunctionBase &f)
adds the function to this card and returns a reference to it.
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > inertia_matrix(const MAST::ElementBase &e) const
std::vector< double > eigen_matrix_to_std_vector(RealMatrixX M)
Converts an Eigen Matrix object to a std::vector.
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_capacitance_matrix(const MAST::ElementBase &e) const
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_expansion_B_matrix(const MAST::ElementBase &e) const
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_expansion_A_matrix(const MAST::ElementBase &e) const
virtual bool depends_on(const MAST::FunctionBase &f) const
returns true if the property card depends on the function f
const MAST::StrainType strain_type() const
returns the type of strain to be used for this element
void set_bending_model(MAST::BendingOperatorType b)
returns the bending model to be used for the 2D element
bool if_diagonal_mass_matrix() const
returns the type of strain to be used for this element
bool plane_stress() const
returns the flag for plane stress.
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_conductance_matrix(const MAST::ElementBase &e) const
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > stiffness_A_matrix(const MAST::ElementBase &e) const
virtual void set_material(MAST::MaterialPropertyCardBase &mat)
sets the material card