5 #include "libmesh/point.h" 21 TEST_CASE(
"element_property_card_constant_heat_transfer_2d",
22 "[heat_transfer],[2D],[isotropic],[constant],[property]")
48 section.
add(thickness_f);
49 section.
add(offset_f);
54 REQUIRE( section.
dim() == dim);
59 SECTION(
"2D section thermal conductance matrix")
73 const libMesh::Point point(2.3, 3.1, 5.2);
74 const Real time = 2.34;
76 conduct_mat->operator()(point, time, D_sec_conduc);
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;
91 CHECK_THAT( test, Catch::Approx<double>(truth) );
96 TEST_CASE(
"element_property_card_constant_transient_heat_transfer_2d",
97 "[heat_transfer],[2D],[isotropic],[constant],[property],[transient]")
126 section.
add(thickness_f);
127 section.
add(offset_f);
132 REQUIRE( section.
dim() == dim);
138 SECTION(
"2D section thermal capacitance matrix")
152 const libMesh::Point point(2.3, 3.1, 5.2);
153 const Real time = 2.34;
155 capaci_mat->operator()(point, time, D_sec_capac);
158 RealMatrixX D_sec_capac_true = RealMatrixX::Zero(1,1);
159 D_sec_capac_true(0,0) = 908.0*1420.5*0.06;
169 CHECK_THAT( test, Catch::Approx<double>(truth) );
174 TEST_CASE(
"element_property_card_constant_thermoelastic_2d",
175 "[thermoelastic],[2D],[isotropic],[constant],[property]")
201 material.
add(alpha_f);
207 section.
add(thickness_f);
208 section.
add(offset_f);
213 REQUIRE( section.
dim() == dim);
218 SECTION(
"2D plane stress thermal expansion A matrix")
239 const libMesh::Point point(2.3, 3.1, 5.2);
240 const Real time = 2.34;
242 texp_A_mat->operator()(point, time, D_sec_texpA);
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;
257 CHECK_THAT( test, Catch::Approx<double>(truth) );
261 SECTION(
"2D plane strain thermal expansion A matrix")
282 const libMesh::Point point(2.3, 3.1, 5.2);
283 const Real time = 2.34;
285 texp_A_mat->operator()(point, time, D_sec_texpA);
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;
300 CHECK_THAT( test, Catch::Approx<double>(truth) );
304 SECTION(
"2D plane stress thermal expansion B matrix")
325 const libMesh::Point point(2.3, 3.1, 5.2);
326 const Real time = 2.34;
328 texp_B_mat->operator()(point, time, D_sec_texpB);
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;
343 CHECK_THAT( test, Catch::Approx<double>(truth) );
346 SECTION(
"2D plane strain thermal expansion B matrix")
367 const libMesh::Point point(2.3, 3.1, 5.2);
368 const Real time = 2.34;
370 texp_B_mat->operator()(point, time, D_sec_texpB);
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;
385 CHECK_THAT( test, Catch::Approx<double>(truth) );
391 "[dynamic],[2D],[isotropic],[constant],[property]")
425 section.
add(thickness_f);
426 section.
add(offset_f);
427 section.
add(kappa_f);
432 REQUIRE( section.
dim() == dim);
444 SECTION(
"2D section inertia matrix")
456 std::unique_ptr<MAST::FieldFunction<RealMatrixX>> inertia_mat = section.
inertia_matrix();
458 const libMesh::Point point(2.3, 3.1, 5.2);
459 const Real time = 2.34;
461 inertia_mat->operator()(point, time, D_sec_iner);
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;
471 D_sec_iner_true *= 1420.5;
481 CHECK_THAT( test, Catch::Approx<double>(truth) );
487 TEST_CASE(
"element_property_card_constant_structural_2d",
488 "[structural],[2D],[isotropic],[constant],[property]")
519 section.
add(thickness_f);
520 section.
add(offset_f);
521 section.
add(kappa_f);
526 REQUIRE( section.
dim() == dim);
533 SECTION(
"solid_2d_section is isotropic")
538 SECTION(
"set_get_strain_type")
552 SECTION(
"set_get_bending_model")
564 SECTION(
"set_get_plane_stress_strain")
586 SECTION(
"2D plane stress extension stiffness matrix")
605 std::unique_ptr<MAST::FieldFunction<RealMatrixX>> extension_stiffness_mat = section.
stiffness_A_matrix();
607 const libMesh::Point point(2.3, 3.1, 5.2);
608 const Real time = 2.34;
610 extension_stiffness_mat->operator()(point, time, D_sec_ext);
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;
628 CHECK_THAT( test, Catch::Approx<double>(truth) );
632 SECTION(
"2D plane strain extension stiffness matrix")
651 std::unique_ptr<MAST::FieldFunction<RealMatrixX>> extension_stiffness_mat = section.
stiffness_A_matrix();
653 const libMesh::Point point(2.3, 3.1, 5.2);
654 const Real time = 2.34;
656 extension_stiffness_mat->operator()(point, time, D_sec_ext);
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;
674 CHECK_THAT( test, Catch::Approx<double>(truth) );
678 SECTION(
"2D plane stress bending section stiffness matrix")
697 std::unique_ptr<MAST::FieldFunction<RealMatrixX>> bending_stiffness_mat = section.
stiffness_D_matrix();
699 const libMesh::Point point(2.3, 3.1, 5.2);
700 const Real time = 2.34;
702 bending_stiffness_mat->operator()(point, time, D_sec_bnd);
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;
720 CHECK_THAT( test, Catch::Approx<double>(truth) );
724 SECTION(
"2D plane strain bending section stiffness matrix")
743 std::unique_ptr<MAST::FieldFunction<RealMatrixX>> bending_stiffness_mat = section.
stiffness_D_matrix();
745 const libMesh::Point point(2.3, 3.1, 5.2);
746 const Real time = 2.34;
748 bending_stiffness_mat->operator()(point, time, D_sec_bnd);
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;
766 CHECK_THAT( test, Catch::Approx<double>(truth) );
770 SECTION(
"2D plane stress extension-bending section stiffness matrix")
790 std::unique_ptr<MAST::FieldFunction<RealMatrixX>> bndext_stiffness_mat = section.
stiffness_B_matrix();
792 const libMesh::Point point(2.3, 3.1, 5.2);
793 const Real time = 2.34;
795 bndext_stiffness_mat->operator()(point, time, D_sec_bndext);
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;
813 CHECK_THAT( test, Catch::Approx<double>(truth) );
817 SECTION(
"2D plane strain extension-bending section stiffness matrix")
837 std::unique_ptr<MAST::FieldFunction<RealMatrixX>> bndext_stiffness_mat = section.
stiffness_B_matrix();
839 const libMesh::Point point(2.3, 3.1, 5.2);
840 const Real time = 2.34;
842 bndext_stiffness_mat->operator()(point, time, D_sec_bndext);
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;
860 CHECK_THAT( test, Catch::Approx<double>(truth) );
864 SECTION(
"2D transverse shear section stiffness matrix")
886 const libMesh::Point point(2.3, 3.1, 5.2);
887 const Real time = 2.34;
889 trans_shear_stiffness_mat->operator()(point, time, D_sec_shr);
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;
904 CHECK_THAT( test, Catch::Approx<double>(truth) );
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...
void set_plane_stress(bool val)
sets the flag for plane stress.
libMesh::LibMeshInit * p_global_init
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