5 #include "libmesh/point.h" 22 TEST_CASE(
"constant_isotropic_thermoelastic_material_1d_sensitivity",
23 "[isotropic],[material],[constant],[1D],[thermoelastic],[sensitivity]")
31 material.
add(alpha_f);
33 REQUIRE( material.
contains(
"alpha_expansion") );
35 SECTION(
"1D material thermal expansion matrix derivative w.r.t. coefficient of thermal expansion")
40 const libMesh::Point point(2.3, 3.1, 5.2);
41 const Real time = 2.34;
43 mat_texp.
derivative(alpha, point, time, dD_texp);
47 dD_texp_true(0,0) = 1.0;
57 CHECK_THAT( test, Catch::Approx<double>(truth) );
62 TEST_CASE(
"constant_isotropic_structural_material_1d_sensitivity",
63 "[isotropic],[material],[constant],[1D],[structural],[sensitivity]")
81 material.
add(kappa_f);
86 REQUIRE( material.
contains(
"kappa") );
88 SECTION(
"1D material constitutive matrix derivative w.r.t. elastic modulus")
93 const libMesh::Point point(2.3, 3.1, 5.2);
94 const Real time = 2.34;
101 dD_true(1,1) = 0.375939849624060;
111 CHECK_THAT( test, Catch::Approx<double>(truth) );
114 SECTION(
"1D material constitutive matrix derivative w.r.t. poisson's ratio")
119 const libMesh::Point point(2.3, 3.1, 5.2);
120 const Real time = 2.34;
122 mat_stiffness.
derivative(nu, point, time, dD);
127 dD_true(1,1) = -0.203516309570920e+11;
137 CHECK_THAT( test, Catch::Approx<double>(truth) );
142 TEST_CASE(
"constant_isotropic_heat_transfer_material_1d_sensitivity",
143 "[isotropic],[material],[1D],[heat_transfer],[constant],[sensitivity]")
158 REQUIRE( material.
contains(
"k_th") );
160 SECTION(
"material depends on the parameters that it should")
165 SECTION(
"material does not depend on other parameters")
171 SECTION(
"1D thermal conductivity matrix derivative w.r.t. thermal conductivity")
176 const libMesh::Point point(2.3, 3.1, 5.2);
177 const Real time = 2.34;
182 RealMatrixX dD_k_true = RealMatrixX::Identity(1,1);
192 CHECK_THAT( test, Catch::Approx<double>(truth) );
197 TEST_CASE(
"constant_isotropic_transient_heat_transfer_material_1d_sensitivity",
198 "[isotropic],[material],[1D],[heat_transfer],[constant],[transient],[sensitivity]")
216 REQUIRE( material.
contains(
"rho") );
219 SECTION(
"material depends on the parameters that it should")
225 SECTION(
"material does not depend on other parameters")
231 SECTION(
"1D thermal capacitance matrix derivative w.r.t. specifc heat capacity")
236 const libMesh::Point point(2.3, 3.1, 5.2);
237 const Real time = 2.34;
239 mat_capacit.
derivative(cp, point, time, dD_cp);
242 RealMatrixX dD_cp_true = RealMatrixX::Identity(1,1);
243 dD_cp_true *= 1234.5;
253 CHECK_THAT( test, Catch::Approx<double>(truth) );
258 TEST_CASE(
"constant_isotropic_thermoelastic_material_2d_sensitivity",
259 "[isotropic],[material],[2D],[thermoelastic],[constant],[sensitivity]")
267 material.
add(alpha_f);
269 REQUIRE( material.
contains(
"alpha_expansion") );
271 SECTION(
"material does not depend on other parameters")
277 SECTION(
"2D material thermal expansion matrix derivative w.r.t. coefficient of thermal expansion")
282 const libMesh::Point point(2.3, 3.1, 5.2);
283 const Real time = 2.34;
285 mat_texp.
derivative(alpha, point, time, dD_texp);
289 dD_texp_true(0,0) = 1.0;
290 dD_texp_true(1,0) = 1.0;
300 CHECK_THAT( test, Catch::Approx<double>(truth) );
305 TEST_CASE(
"constant_isotropic_structural_material_2d_sensitivity",
306 "[isotropic],[material],[2D],[structural],[constant],[sensitivity]")
324 material.
add(kappa_f);
329 REQUIRE( material.
contains(
"kappa") );
331 SECTION(
"material depends on the parameters that it should")
338 SECTION(
"material does not depend on other parameters")
344 SECTION(
"2D plane stress material constitutive matrix derivative w.r.t. elastic modulus is correct")
346 const bool is_plane_stress =
true;
350 const libMesh::Point point(2.3, 3.1, 5.2);
351 const Real time = 2.34;
353 mat_stiffness.
derivative(E, point, time, dK_mat);
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;
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());
371 CHECK_THAT( test, Catch::Approx<double>(truth) );
375 SECTION(
"2D plane stress material constitutive matrix derivative w.r.t. poisson's ratio is correct")
377 const bool is_plane_stress =
true;
381 const libMesh::Point point(2.3, 3.1, 5.2);
382 const Real time = 2.34;
384 mat_stiffness.
derivative(nu, point, time, dK_mat);
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;
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());
402 CHECK_THAT( test, Catch::Approx<double>(truth) );
405 SECTION(
"2D plane strain material constitutive matrix derivative w.r.t. elastic modulus is correct")
407 const bool is_plane_stress =
false;
411 const libMesh::Point point(2.3, 3.1, 5.2);
412 const Real time = 2.34;
414 mat_stiffness.
derivative(E, point, time, dK_mat);
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;
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());
453 CHECK_THAT( test, Catch::Approx<double>(truth) );
456 SECTION(
"2D plane strain material constitutive matrix derivative w.r.t. poisson's ratio is correct")
458 const bool is_plane_stress =
false;
462 const libMesh::Point point(2.3, 3.1, 5.2);
463 const Real time = 2.34;
465 mat_stiffness.
derivative(nu, point, time, dK_mat);
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;
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());
483 CHECK_THAT( test, Catch::Approx<double>(truth) );
486 SECTION(
"transverse shear material stiffness matrix derivative w.r.t. elastic modulus")
491 const libMesh::Point point(2.3, 3.1, 5.2);
492 const Real time = 2.34;
494 mat_trans_shear.
derivative(E, point, time, dD_trans_shear);
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;
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());
509 CHECK_THAT( test, Catch::Approx<double>(truth) );
512 SECTION(
"transverse shear material stiffness matrix derivative w.r.t. poisson's ratio")
517 const libMesh::Point point(2.3, 3.1, 5.2);
518 const Real time = 2.34;
520 mat_trans_shear.
derivative(nu, point, time, dD_trans_shear);
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;
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());
535 CHECK_THAT( test, Catch::Approx<double>(truth) );
538 SECTION(
"transverse shear material stiffness matrix derivative w.r.t. shear coefficient")
543 const libMesh::Point point(2.3, 3.1, 5.2);
544 const Real time = 2.34;
546 mat_trans_shear.
derivative(kappa, point, time, dD_trans_shear);
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;
561 CHECK_THAT( test, Catch::Approx<double>(truth) );
566 TEST_CASE(
"constant_isotropic_transient_heat_transfer_material_2d_sensitivity",
567 "[isotropic],[material],[2D],[heat_transfer],[constant],[sensitivity],[transient]")
588 REQUIRE( material.
contains(
"rho") );
590 REQUIRE( material.
contains(
"k_th") );
592 SECTION(
"material depends on the parameters that it should")
599 SECTION(
"material does not depend on other parameters")
605 SECTION(
"2D thermal capacitance matrix derivative w.r.t. specifc heat capacity")
610 const libMesh::Point point(2.3, 3.1, 5.2);
611 const Real time = 2.34;
613 mat_capacit.
derivative(cp, point, time, dD_cp);
616 RealMatrixX dD_cp_true = RealMatrixX::Identity(1,1);
617 dD_cp_true *= 1234.5;
627 CHECK_THAT( test, Catch::Approx<double>(truth) );
630 SECTION(
"2D thermal capacitance matrix derivative w.r.t. density")
635 const libMesh::Point point(2.3, 3.1, 5.2);
636 const Real time = 2.34;
638 mat_capacit.
derivative(rho, point, time, dD_cp);
641 RealMatrixX dD_cp_true = RealMatrixX::Identity(1,1);
652 CHECK_THAT( test, Catch::Approx<double>(truth) );
657 TEST_CASE(
"constant_isotropic_heat_transfer_material_2d_sensitivity",
658 "[isotropic],[material],[2D],[heat_transfer][constant],[sensitivity]")
679 REQUIRE( material.
contains(
"rho") );
681 REQUIRE( material.
contains(
"k_th") );
683 SECTION(
"material depends on the parameters that it should")
690 SECTION(
"material does not depend on other parameters")
696 SECTION(
"2D thermal conductivity matrix derivative w.r.t. thermal conductivity")
701 const libMesh::Point point(2.3, 3.1, 5.2);
702 const Real time = 2.34;
707 RealMatrixX dD_k_true = RealMatrixX::Identity(2,2);
717 CHECK_THAT( test, Catch::Approx<double>(truth) );
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...
virtual void derivative(const MAST::FunctionBase &f, ValType &v) const
calculates the value of the function derivative and returns it in v.
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