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