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