5 #include "libmesh/point.h" 18 TEST_CASE(
"constant_isotropic_thermoelastic_material_1d",
19 "[isotropic],[material],[constant],[1D],[thermoelastic]")
27 material.
add(alpha_f);
29 REQUIRE( material.
contains(
"alpha_expansion") );
31 SECTION(
"1D material thermal expansion matrix")
36 const libMesh::Point point(2.3, 3.1, 5.2);
37 const Real time = 2.34;
39 mat_texp(point, time, D_texp);
43 D_texp_true(0,0) = 5.43e-05;
53 CHECK_THAT( test, Catch::Approx<double>(truth) );
58 TEST_CASE(
"constant_isotropic_structural_material_1d",
59 "[isotropic],[material],[constant],[1D],[structural]")
80 SECTION(
"1D material constitutive matrix")
85 const libMesh::Point point(2.3, 3.1, 5.2);
86 const Real time = 2.34;
88 mat_stiffness(point, time, D);
93 D_true(1,1) = 2.706766917293233e+10;
103 CHECK_THAT( test, Catch::Approx<double>(truth) );
108 TEST_CASE(
"constant_isotropic_heat_transfer_material_1d",
109 "[isotropic],[material],[1D],[heat_transfer][constant]")
124 REQUIRE( material.
contains(
"k_th") );
126 SECTION(
"material depends on the parameters that it should")
131 SECTION(
"material does not depend on other parameters")
137 SECTION(
"1D thermal conductivity matrix")
142 const libMesh::Point point(2.3, 3.1, 5.2);
143 const Real time = 2.34;
145 mat_conduct(point, time, D_k);
159 CHECK_THAT( test, Catch::Approx<double>(truth) );
164 TEST_CASE(
"constant_isotropic_transient_heat_transfer_material_1d",
165 "[isotropic],[material],[1D],[heat_transfer],[constant],[transient]")
183 REQUIRE( material.
contains(
"rho") );
186 SECTION(
"material depends on the parameters that it should")
192 SECTION(
"material does not depend on other parameters")
198 SECTION(
"1D thermal capacitance matrix")
203 const libMesh::Point point(2.3, 3.1, 5.2);
204 const Real time = 2.34;
206 mat_capacit(point, time, D_cp);
209 RealMatrixX D_cp_true = RealMatrixX::Identity(1,1);
210 D_cp_true *= (908.0 * 1234.5);
220 CHECK_THAT( test, Catch::Approx<double>(truth) );
225 TEST_CASE(
"constant_isotropic_thermoelastic_material_2d",
226 "[isotropic],[material],[2D],[thermoelastic][constant]")
234 material.
add(alpha_f);
236 REQUIRE( material.
contains(
"alpha_expansion") );
238 SECTION(
"material does not depend on other parameters")
244 SECTION(
"2D material thermal expansion matrix")
249 const libMesh::Point point(2.3, 3.1, 5.2);
250 const Real time = 2.34;
252 mat_texp(point, time, D_texp);
256 D_texp_true(0,0) = 5.43e-05;
257 D_texp_true(1,0) = 5.43e-05;
267 CHECK_THAT( test, Catch::Approx<double>(truth) );
273 "[isotropic],[material],[2D],[structural],[constant]")
279 const Real G = E() / (2.0 * (1.0+nu()));
296 SECTION(
"material depends on the parameters that it should")
302 SECTION(
"material does not depend on other parameters")
308 SECTION(
"2D plane stress material constitutive matrix is correct")
310 const bool is_plane_stress =
true;
314 const libMesh::Point point(2.3, 3.1, 5.2);
315 const Real time = 2.34;
317 mat_stiffness(point, time, K_mat);
329 REQUIRE_THAT( test, Catch::Approx<double>(truth) );
332 SelfAdjointEigenSolver<RealMatrixX> eigensolver(K_mat);
333 if (eigensolver.info() != Success)
335 std::cout <<
"eigensolver failed to converge!" << std::endl;
338 RealVectorX eigenvalues = eigensolver.eigenvalues();
339 for (uint i=0; i<K_mat.cols(); i++)
341 REQUIRE( eigenvalues(i) > 0.0 );
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;
360 CHECK_THAT( test, Catch::Approx<double>(truth) );
363 SECTION(
"2D plane strain material constitutive matrix is correct")
365 const bool is_plane_stress =
false;
369 const libMesh::Point point(2.3, 3.1, 5.2);
370 const Real time = 2.34;
372 mat_stiffness(point, time, K_mat);
384 REQUIRE_THAT( test, Catch::Approx<double>(truth) );
387 SelfAdjointEigenSolver<RealMatrixX> eigensolver(K_mat);
388 if (eigensolver.info() != Success)
390 std::cout <<
"eigensolver failed to converge!" << std::endl;
393 RealVectorX eigenvalues = eigensolver.eigenvalues();
394 for (uint i=0; i<K_mat.cols(); i++)
396 REQUIRE( eigenvalues(i) > 0.0 );
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;
415 CHECK_THAT( test, Catch::Approx<double>(truth) );
418 SECTION(
"transverse shear material stiffness matrix")
423 const libMesh::Point point(2.3, 3.1, 5.2);
424 const Real time = 2.34;
426 mat_trans_shear(point, time, D_trans_shear);
429 RealMatrixX D_trans_shear_true = RealMatrixX::Zero(2,2);
430 D_trans_shear_true(0,0) = G;
431 D_trans_shear_true(1,1) = G;
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());
441 CHECK_THAT( test, Catch::Approx<double>(truth) );
446 TEST_CASE(
"constant_isotropic_heat_transfer_material_2d",
447 "[isotropic],[material],[2D],[heat_transfer][constant]")
462 REQUIRE( material.
contains(
"k_th") );
464 SECTION(
"material depends on the parameters that it should")
469 SECTION(
"material does not depend on other parameters")
475 SECTION(
"2D thermal conductivity matrix")
480 const libMesh::Point point(2.3, 3.1, 5.2);
481 const Real time = 2.34;
483 mat_conduct(point, time, D_k);
497 CHECK_THAT( test, Catch::Approx<double>(truth) );
502 TEST_CASE(
"constant_isotropic_transient_heat_transfer_material_2d",
503 "[isotropic],[material],[2D],[heat_transfer],[constant],[transient]")
521 REQUIRE( material.
contains(
"rho") );
524 SECTION(
"material depends on the parameters that it should")
530 SECTION(
"material does not depend on other parameters")
536 SECTION(
"2D thermal capacitance matrix")
541 const libMesh::Point point(2.3, 3.1, 5.2);
542 const Real time = 2.34;
544 mat_capacit(point, time, D_cp);
547 RealMatrixX D_cp_true = RealMatrixX::Identity(1,1);
548 D_cp_true *= (908.0 * 1234.5);
558 CHECK_THAT( test, Catch::Approx<double>(truth) );
563 TEST_CASE(
"constant_isotropic_thermoelastic_material_3d",
564 "[isotropic],[material],[3D],[thermoelastic][constant]")
572 material.
add(alpha_f);
574 REQUIRE( material.
contains(
"alpha_expansion") );
576 SECTION(
"material does not depend on other parameters")
582 SECTION(
"3D material thermal expansion matrix")
587 const libMesh::Point point(2.3, 3.1, 5.2);
588 const Real time = 2.34;
590 mat_texp(point, time, D_texp);
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;
606 CHECK_THAT( test, Catch::Approx<double>(truth) );
612 "[isotropic],[material],[3D],[structural],[constant]")
618 const Real G = E() / (2.0 * (1.0+nu()));
635 SECTION(
"material depends on the parameters that it should")
641 SECTION(
"material does not depend on other parameters")
647 SECTION(
"3D material constitutive matrix is correct")
649 const bool is_plane_stress =
true;
653 const libMesh::Point point(2.3, 3.1, 5.2);
654 const Real time = 2.34;
656 mat_stiffness(point, time, K_mat);
668 REQUIRE_THAT( test, Catch::Approx<double>(truth) );
671 SelfAdjointEigenSolver<RealMatrixX> eigensolver(K_mat);
672 if (eigensolver.info() != Success)
674 std::cout <<
"eigensolver failed to converge!" << std::endl;
677 RealVectorX eigenvalues = eigensolver.eigenvalues();
678 for (uint i=0; i<K_mat.cols(); i++)
680 REQUIRE( eigenvalues(i) > 0.0 );
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.;
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;
705 CHECK_THAT( test, Catch::Approx<double>(truth) );
710 TEST_CASE(
"constant_isotropic_heat_transfer_material_3d",
711 "[isotropic],[material],[3D],[heat_transfer][constant]")
726 REQUIRE( material.
contains(
"k_th") );
728 SECTION(
"material depends on the parameters that it should")
733 SECTION(
"material does not depend on other parameters")
739 SECTION(
"2D thermal conductivity matrix")
744 const libMesh::Point point(2.3, 3.1, 5.2);
745 const Real time = 2.34;
747 mat_conduct(point, time, D_k);
761 CHECK_THAT( test, Catch::Approx<double>(truth) );
766 TEST_CASE(
"constant_isotropic_transient_heat_transfer_material_3d",
767 "[isotropic],[material],[3D],[heat_transfer],[constant],[transient]")
785 REQUIRE( material.
contains(
"rho") );
788 SECTION(
"material depends on the parameters that it should")
794 SECTION(
"material does not depend on other parameters")
800 SECTION(
"2D thermal capacitance matrix")
805 const libMesh::Point point(2.3, 3.1, 5.2);
806 const Real time = 2.34;
808 mat_capacit(point, time, D_cp);
811 RealMatrixX D_cp_true = RealMatrixX::Identity(1,1);
812 D_cp_true *= (908.0 * 1234.5);
822 CHECK_THAT( test, Catch::Approx<double>(truth) );
828 "[isotropic],[material],[3D],[dynamic],[constant]")
843 REQUIRE( material.
contains(
"rho") );
845 SECTION(
"material depends on the parameters that it should")
850 SECTION(
"material does not depend on other parameters")
856 SECTION(
"3D inertia matrix")
861 const libMesh::Point point(2.3, 3.1, 5.2);
862 const Real time = 2.34;
864 mat_capacit(point, time, D_inertia);
867 RealMatrixX D_inertia_true = RealMatrixX::Identity(3,3);
868 D_inertia_true *= 1234.5;
878 CHECK_THAT( test, Catch::Approx<double>(truth) );
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...
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