MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_quad4_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("quad4_linear_structural_inertial_consistent",
46  "[quad],[quad4],[dynamic],[2D],[element]")
47 {
48  RealMatrixX coords = RealMatrixX::Zero(3,4);
49  coords << -1.0, 1.0, 1.0, -1.0,
50  -1.0, -1.0, 1.0, 1.0,
51  0.0, 0.0, 0.0, 0.0;
52  TEST::TestStructuralSingleElement2D test_elem(libMesh::QUAD4, coords);
53 
54  const Real V0 = test_elem.reference_elem->volume();
55  REQUIRE(test_elem.reference_elem->volume() == TEST::get_shoelace_area(coords));
56 
57  // Set the strain type to linear for the section
59 
60  // Set the bending operator
62 
63  // Set mass matrix to consistent & compute inertial residual and Jacobians.
64  test_elem.section.set_diagonal_mass_matrix(false);
65  REQUIRE_FALSE(test_elem.section.if_diagonal_mass_matrix());
67 
68  double val_margin = (test_elem.jacobian_xddot0.array().abs()).mean() * 1.490116119384766e-08;
69 
70  libMesh::out << "Jac_xddot0:\n" << test_elem.jacobian_xddot0 << std::endl;
71 
72 
73  SECTION("Inertial Jacobian (mass matrix) finite difference check")
74  {
76  test_elem.elem_solution, test_elem.jacobian_fd);
77 
78  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
79 
80  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0),
81  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
82  }
83 
84 
85  SECTION("Inertial Jacobian (mass matrix) should be symmetric")
86  {
87  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0),
88  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0.transpose())));
89  }
90 
91 
92  SECTION("Determinant of undeformed inertial Jacobian (mass matrix) should be zero")
93  {
94  REQUIRE(test_elem.jacobian_xddot0.determinant() == Approx(0.0).margin(1e-06));
95  }
96 
97 
98  SECTION("Inertial Jacobian (mass matrix) eigenvalues should all be positive")
99  {
100  SelfAdjointEigenSolver<RealMatrixX> eigensolver(test_elem.jacobian_xddot0, false);
101  RealVectorX eigenvalues = eigensolver.eigenvalues();
102  libMesh::out << "Eigenvalues are:\n" << eigenvalues << std::endl;
103  REQUIRE(eigenvalues.minCoeff()>0.0);
104  }
105 }
106 
107 
108 TEST_CASE("quad4_linear_structural_inertial_lumped",
109  "[quad],[quad4],[dynamic],[2D],[element]")
110 {
111  RealMatrixX coords = RealMatrixX::Zero(3,4);
112  coords << -1.0, 1.0, 1.0, -1.0,
113  -1.0, -1.0, 1.0, 1.0,
114  0.0, 0.0, 0.0, 0.0;
115  TEST::TestStructuralSingleElement2D test_elem(libMesh::QUAD4, coords);
116 
117  const Real V0 = test_elem.reference_elem->volume();
118  REQUIRE(test_elem.reference_elem->volume() == TEST::get_shoelace_area(coords));
119 
120  // Set the strain type to linear for the section
122 
123  // Set the bending operator
125 
126  // Set mass matrix to lumped & compute inertial residual and Jacobians.
127  test_elem.section.set_diagonal_mass_matrix(true);
128  REQUIRE(test_elem.section.if_diagonal_mass_matrix());
130 
131  double val_margin = (test_elem.jacobian_xddot0.array().abs()).mean() * 1.490116119384766e-08;
132 
133  libMesh::out << "Jac_xddot0:\n" << test_elem.jacobian_xddot0 << std::endl;
134 
135 
136  SECTION("Inertial Jacobian (mass matrix) finite difference check")
137  {
139  test_elem.elem_solution, test_elem.jacobian_fd);
140 
141  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
142 
143  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0),
144  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
145  }
146 
147 
148  SECTION("Inertial Jacobian (mass matrix) should be symmetric")
149  {
150  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0),
151  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_xddot0.transpose())));
152  }
153 
154 
155  SECTION("Determinant of undeformed inertial Jacobian (mass matrix) should be zero")
156  {
157  REQUIRE(test_elem.jacobian_xddot0.determinant() == Approx(0.0).margin(1e-06));
158  }
159 
160 
161  SECTION("Inertial Jacobian (mass matrix) eigenvalues should all be positive")
162  {
163  SelfAdjointEigenSolver<RealMatrixX> eigensolver(test_elem.jacobian_xddot0, false);
164  RealVectorX eigenvalues = eigensolver.eigenvalues();
165  libMesh::out << "Eigenvalues are:\n" << eigenvalues << std::endl;
166  REQUIRE(eigenvalues.minCoeff()>0.0);
167  }
168 }
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...
libMesh::Elem * reference_elem
Pointer to the actual libMesh element object.
Definition: mast_mesh.h:51
RealMatrixX jacobian_xddot0
Matrix storage for acceleration Jacobian of baseline/undeformed element.
void set_diagonal_mass_matrix(bool m)
sets the mass matrix to be diagonal or consistent
Real get_shoelace_area(RealMatrixX X)
Calcualtes the area of a 2D polygon using the shoelace formula.
MAST::Solid2DSectionElementPropertyCard section
libMesh::Real Real
void set_strain(MAST::StrainType strain)
sets the type of strain to be used, which is LINEAR_STRAIN by default
libMesh::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26
std::vector< double > eigen_matrix_to_std_vector(RealMatrixX M)
Converts an Eigen Matrix object to a std::vector.
Matrix< Real, Dynamic, Dynamic > RealMatrixX
TEST_CASE("quad4_linear_structural_inertial_consistent", "[quad],[quad4],[dynamic],[2D],[element]")
Matrix< Real, Dynamic, 1 > RealVectorX
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
RealMatrixX jacobian_fd
Matrix storage for element Jacobian approximated by finite difference.