MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_edge2_linear_structural_inertial_jacobian.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2020 Manav Bhatia and MAST authors
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 // libMesh includes
21 #include "libmesh/libmesh.h"
22 #include "libmesh/elem.h"
23 #include "libmesh/dof_map.h"
24 
25 // MAST includes
26 #include "base/parameter.h"
34 #include "base/nonlinear_system.h"
35 #include "mesh/geom_elem.h"
36 
37 // Test includes
38 #include "catch.hpp"
39 #include "test_helpers.h"
41 
42 extern libMesh::LibMeshInit* p_global_init;
43 
44 
45 TEST_CASE("edge2_linear_structural_inertial_consistent",
46  "[1D],[dynamic],[edge],[edge2],[element]")
47 {
48  RealMatrixX coords = RealMatrixX::Zero(3, 2);
49  coords << -1.0, 1.0, 0.0,
50  0.0, 0.0, 0.0;
51  TEST::TestStructuralSingleElement1D test_elem(libMesh::EDGE2, coords);
52 
53  const Real V0 = test_elem.reference_elem->volume();
54 
55  // Set mass matrix to consistent & compute inertial residual and Jacobians.
56  test_elem.section.set_diagonal_mass_matrix(false);
58 
59  double val_margin = (test_elem.jacobian_xddot0.array().abs()).mean() * 1.490116119384766e-08;
60 
61  libMesh::out << "Jac_xddot0:\n" << test_elem.jacobian_xddot0 << std::endl;
62 
63 
64  SECTION("Inertial Jacobian (mass matrix) finite difference check")
65  {
67  test_elem.elem_solution, test_elem.jacobian_fd);
68 
69  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
70 
71  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0),
72  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
73  }
74 
75 
76  SECTION("Inertial Jacobian (mass matrix) should be symmetric")
77  {
78  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0),
79  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0.transpose())));
80  }
81 
82 
83  SECTION("Determinant of undeformed inertial Jacobian (mass matrix) should be zero")
84  {
85  REQUIRE(test_elem.jacobian_xddot0.determinant() == Approx(0.0).margin(1e-06));
86  }
87 
88 
89  SECTION("Inertial Jacobian (mass matrix) eigenvalues should all be positive")
90  {
91  SelfAdjointEigenSolver<RealMatrixX> eigensolver(test_elem.jacobian_xddot0, false);
92  RealVectorX eigenvalues = eigensolver.eigenvalues();
93  //libMesh::out << "Eigenvalues are:\n" << eigenvalues << std::endl;
94  REQUIRE(eigenvalues.minCoeff()>0.0);
95  }
96 }
97 
98 
99 TEST_CASE("edge2_linear_structural_inertial_lumped",
100  "[1D],[dynamic],[edge],[edge2]")
101 {
102  RealMatrixX coords = RealMatrixX::Zero(3, 2);
103  coords << -1.0, 1.0, 0.0,
104  0.0, 0.0, 0.0;
105  TEST::TestStructuralSingleElement1D test_elem(libMesh::EDGE2, coords);
106 
107  const Real V0 = test_elem.reference_elem->volume();
108 
109  // Set mass matrix to lumped & compute inertial residual and Jacobians.
110  test_elem.section.set_diagonal_mass_matrix(true);
112 
113  double val_margin = (test_elem.jacobian_xddot0.array().abs()).mean() * 1.490116119384766e-08;
114 
115  libMesh::out << "Jac_xddot0:\n" << test_elem.jacobian_xddot0 << std::endl;
116 
117 
118  SECTION("Inertial Jacobian (mass matrix) finite difference check")
119  {
121  test_elem.elem_solution, test_elem.jacobian_fd);
122 
123  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
124 
125  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0),
126  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
127  }
128 
129 
130  SECTION("Inertial Jacobian (mass matrix) should be symmetric")
131  {
132  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0),
133  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0.transpose())));
134  }
135 
136 
137  SECTION("Determinant of undeformed inertial Jacobian (mass matrix) should be zero")
138  {
139  REQUIRE(test_elem.jacobian_xddot0.determinant() == Approx(0.0).margin(1e-06));
140  }
141 
142 
143  SECTION("Inertial Jacobian (mass matrix) eigenvalues should all be positive")
144  {
145  SelfAdjointEigenSolver<RealMatrixX> eigensolver(test_elem.jacobian_xddot0, false);
146  RealVectorX eigenvalues = eigensolver.eigenvalues();
147  //libMesh::out << "Eigenvalues are:\n" << eigenvalues << std::endl;
148  REQUIRE(eigenvalues.minCoeff()>0.0);
149  }
150 }
void approximate_inertial_jacobian_with_finite_difference(MAST::StructuralElementBase &elem, const RealVectorX &initial_elem_solution, RealMatrixX &jacobian)
Approximates the inertial jacobian using a 6th order accurate central finite difference scheme...
TEST_CASE("edge2_linear_structural_inertial_consistent", "[1D],[dynamic],[edge],[edge2],[element]")
RealMatrixX jacobian_xddot0
Matrix storage for acceleration Jacobian of baseline/undeformed element.
libMesh::Elem * reference_elem
Pointer to the actual libMesh element object.
Definition: mast_mesh.h:51
void set_diagonal_mass_matrix(bool m)
sets the mass matrix to be diagonal or consistent
MAST::Solid1DSectionElementPropertyCard section
libMesh::Real Real
std::vector< double > eigen_matrix_to_std_vector(RealMatrixX M)
Converts an Eigen Matrix object to a std::vector.
Matrix< Real, Dynamic, Dynamic > RealMatrixX
RealMatrixX jacobian_fd
Matrix storage for element Jacobian approximated by finite difference.
Matrix< Real, Dynamic, 1 > RealVectorX
libMesh::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26