MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_edge2_linear_structural_internal_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",
46  "[1D],[structural],[edge],[edge2],[linear]")
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 
56 
57  double val_margin = (test_elem.jacobian0.array().abs()).mean() * 1.490116119384766e-08;
58 
59  libMesh::out << "J =\n" << test_elem.jacobian0 << std::endl;
60 
61 
62  SECTION("Internal Jacobian (stiffness matrix) finite difference check")
63  {
65  test_elem.elem_solution, test_elem.jacobian_fd);
66 
67  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
68 
69  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
70  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
71  }
72 
73 
74  SECTION("Internal Jacobian (stiffness matrix) should be symmetric")
75  {
76  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
77  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0.transpose())));
78  }
79 
80 
81  SECTION("Determinant of undeformed internal Jacobian (stiffness matrix) should be zero")
82  {
83  REQUIRE(test_elem.jacobian0.determinant() == Approx(0.0).margin(1e-06));
84  }
85 
86 
87  SECTION("Internal Jacobian (stiffness matrix) eigenvalue check")
88  {
89  /*
90  * Number of zero eigenvalues should equal the number of rigid body
91  * modes. With extension (including torsion) and bending about y and
92  * z axes, we have 6 rigid body modes (3 translations, 3 rotations).
93  *
94  * Note that the use of reduced integration can result in more rigid
95  * body modes than expected.
96  */
97  SelfAdjointEigenSolver<RealMatrixX> eigensolver(test_elem.jacobian0, false);
98  RealVectorX eigenvalues = eigensolver.eigenvalues();
99  libMesh::out << "Eigenvalues are:\n" << eigenvalues << std::endl;
100  uint nz = 0;
101  for (uint i=0; i<eigenvalues.size(); i++)
102  {
103  if (std::abs(eigenvalues(i))<0.0001220703125)
104  {
105  nz++;
106  }
107  }
108  REQUIRE( nz == 6);
109 
113  REQUIRE(eigenvalues.minCoeff()>(-0.0001220703125));
114  }
115 
116 
117  SECTION("Internal Jacobian (stiffness matrix) is invariant to section orientation")
118  {
119  test_elem.section.clear();
120  RealVectorX orientation = RealVectorX::Zero(3);
121  orientation(2) = 1.0;
122  test_elem.section.y_vector() = orientation;
123  test_elem.section.init();
124  test_elem.update_residual_and_jacobian();
125 
126  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
127  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
128  }
129 
130 
131  SECTION("Internal Jacobian (stiffness matrix) is invariant to displacement solution")
132  {
133  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
134  elem_sol << 0.05727841, 0.08896581, 0.09541619, -0.03774913,
135  0.07510557, -0.07122266, -0.00979117, -0.08300009,
136  -0.03453369, -0.05487761, -0.01407677, -0.09268421;
137  test_elem.elem->set_solution(elem_sol);
138  test_elem.update_residual_and_jacobian();
139 
140  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
141  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
142  }
143 
144 
145  SECTION("Internal Jacobian (stiffness matrix) is invariant to element x-location")
146  {
147  TEST::transform_element(test_elem.mesh, coords,5.2, 0.0, 0.0,
148  1.0, 1.0, 0.0, 0.0, 0.0);
149  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
150  test_elem.update_residual_and_jacobian();
151 
152  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
153  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
154  }
155 
156 
157  SECTION("Internal Jacobian (stiffness matrix) is invariant to element y-location")
158  {
159  TEST::transform_element(test_elem.mesh, coords, 0.0, -11.5, 0.0,
160  1.0, 1.0, 0.0, 0.0, 0.0);
161  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
162  test_elem.update_residual_and_jacobian();
163 
164  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
165  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
166  }
167 
168 
169  SECTION("Internal Jacobian (stiffness matrix) is invariant to element z-location")
170  {
171  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 7.6,
172  1.0, 1.0, 0.0, 0.0, 0.0);
173  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
174  test_elem.update_residual_and_jacobian();
175 
176  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
177  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
178  }
179 
180 
181  SECTION("Internal Jacobian (stiffness matrix) checks for element aligned along y-axis")
182  {
183  /*
184  * NOTE: We could try to use the transform_element method here, but the
185  * issue is that if the sin and cos calculations are not exact, then we
186  * may not be perfectly aligned along the y axis like we want.
187  */
188  RealMatrixX new_coordinates = RealMatrixX::Zero(3,test_elem.n_nodes);
189  new_coordinates << 0.0, 0.0, -1.0, 1.0, 0.0, 0.0;
190  test_elem.update_coordinates(new_coordinates);
191  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
192  test_elem.update_residual_and_jacobian();
193 
194  // Finite difference Jacobian check
196  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
197  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
198  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
199 
200  // Symmetry check
201  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
202  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
203 
204  // Determinant check
205  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
206  }
207 
208 
209  SECTION("Internal Jacobian (stiffness matrix) checks for element aligned along z-axis")
210  {
211  /*
212  * NOTE: We could try to use the transform_element method here, but the
213  * issue is that if the sin and cos calculations are not exact, then we
214  * may not be perfectly aligned along the z axis like we want.
215  */
216  RealMatrixX new_coordinates = RealMatrixX::Zero(3,test_elem.n_nodes);
217  new_coordinates << 0.0, 0.0, 0.0, 0.0, -1.0, 1.0;
218  test_elem.update_coordinates(new_coordinates);
219  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
220  test_elem.update_residual_and_jacobian();
221 
222 
223  // Finite difference Jacobian check
225  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
226  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
227  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
228 
229  // Symmetry check
230  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
231  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
232 
233  // Determinant check
234  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
235  }
236 
237 
238  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about z-axis")
239  {
240  // Rotated 63.4 about z-axis at element's centroid
241  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
242  1.0, 1.0, 0.0, 0.0, 63.4);
243  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
244  test_elem.update_residual_and_jacobian();
245 
246  // Finite difference Jacobian check
248  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
249  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
250  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
251 
252  // Symmetry check
253  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
254  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
255 
256  // Determinant check
257  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
258  }
259 
260 
261  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about y-axis")
262  {
263  // Rotated 35.8 about y-axis at element's centroid
264  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
265  1.0, 1.0, 0.0, 35.8, 0.0);
266  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
267  test_elem.update_residual_and_jacobian();
268 
269  // Finite difference Jacobian check
271  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
272  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
273  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
274 
275  // Symmetry check
276  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
277  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
278 
279  // Determinant check
280  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
281  }
282 
283 
284  SECTION("Internal Jacobian (stiffness matrix) checks for element stretched in x-direction")
285  {
286  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
287  3.2, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0);
288  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
289  test_elem.update_residual_and_jacobian();
290 
291  // Finite difference Jacobian check
293  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
294  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
295  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
296 
297  // Symmetry check
298  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
299  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
300 
301  // Determinant check
302  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
303  }
304 
305 
306  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, and rotated")
307  {
308  // Apply arbitrary transformation to the element
309  TEST::transform_element(test_elem.mesh, coords, -5.0, 7.8, -13.1,
310  2.7, 6.4, 20.0, 47.8, -70.1);
311  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
312  test_elem.update_residual_and_jacobian();
313 
314  // Finite difference Jacobian check
316  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
317  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
318  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
319 
320  // Symmetry check
321  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
322  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
323 
324  // Determinant check
325  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
326  }
327 
328 
329  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and displaced")
330  {
331  // Apply arbitrary transformation to the element
332  TEST::transform_element(test_elem.mesh, coords, 4.1, -6.3, 7.5,
333  4.2, 1.5, -18.0, -24.8, 30.1);
334 
335  // Apply arbitrary displacement to the element
336  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
337  elem_sol << 0.08158724, 0.07991906, -0.00719128, 0.02025461,
338  -0.04602193, 0.05280159, 0.03700081, 0.04636344,
339  0.05559377, 0.06448206, 0.08919238, -0.03079122;
340  test_elem.elem->set_solution(elem_sol);
341 
342  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
343  test_elem.update_residual_and_jacobian();
344 
345  // Finite difference Jacobian check
347  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
348  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
349  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
350 
351  // Symmetry check
352  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
353  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
354 
355  // Determinant check
356  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
357  }
358 }
RealMatrixX jacobian
Matrix storage for Jacobian of the element in a perturbed/modified state.
libMesh::Elem * reference_elem
Pointer to the actual libMesh element object.
Definition: mast_mesh.h:51
RealVectorX & y_vector()
returns value of the property val.
MAST::Solid1DSectionElementPropertyCard section
libMesh::Real Real
TEST_CASE("edge2_linear_structural", "[1D],[structural],[edge],[edge2],[linear]")
int n_nodes
Number of nodes per element in the test mesh.
Definition: mast_mesh.h:49
void transform_element(libMesh::MeshBase &mesh, const RealMatrixX X0, Real shift_x, Real shift_y, Real shift_z, Real scale_x, Real scale_y, Real rotation_x, Real rotation_y, Real rotation_z, Real shear_x=0, Real shear_y=0)
Transform an element by applying any combination of: shifts, scales, rotations, and shears...
std::vector< double > eigen_matrix_to_std_vector(RealMatrixX M)
Converts an Eigen Matrix object to a std::vector.
void approximate_internal_jacobian_with_finite_difference(MAST::StructuralElementBase &elem, const RealVectorX &initial_elem_solution, RealMatrixX &jacobian)
Approximates the internal Jacobian of an element using a 6th order accurate central finite difference...
RealMatrixX jacobian0
Matrix storage for Jacobian of baseline/undeformed element.
Matrix< Real, Dynamic, Dynamic > RealMatrixX
RealMatrixX jacobian_fd
Matrix storage for element Jacobian approximated by finite difference.
virtual void set_solution(const RealVectorX &vec, bool if_sens=false)
stores vec as solution for element level calculations, or its sensitivity if if_sens is true...
Matrix< Real, Dynamic, 1 > RealVectorX
libMesh::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26
void update_coordinates(RealMatrixX &new_coordinates)
Update the nodal coordinates in the mesh.
Definition: mast_mesh.h:113
libMesh::ReplicatedMesh mesh
The actual libMesh mesh object.
Definition: mast_mesh.h:52