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