MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_edge2_linear_structural_extension_bending_coupling_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_coupling_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_yy = 0.0;
57  test_elem.kappa_zz = 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  Real eps = 1.4901161193847656e-07;
103  SelfAdjointEigenSolver<RealMatrixX> eigensolver(test_elem.jacobian0, false);
104  RealVectorX eigenvalues = eigensolver.eigenvalues();
105  libMesh::out << "Eigenvalues are:\n" << eigenvalues << std::endl;
106  uint nz = 0;
107  for (uint i=0; i<eigenvalues.size(); i++)
108  {
109  if (std::abs(eigenvalues(i))<eps)
110  {
111  nz++;
112  }
113  }
114  libMesh::out << "Eigenvalues:\n" << eigenvalues << std::endl;
115  REQUIRE( nz == 6);
116 
117  /*
118  * All non-zero eigenvalues should be positive.
119  */
120  REQUIRE(eigenvalues.minCoeff()>(-eps));
121  }
122 
123 
124  SECTION("Internal Jacobian (stiffness matrix) is invariant to section orientation")
125  {
126  test_elem.section.clear();
127  RealVectorX orientation = RealVectorX::Zero(3);
128  orientation(2) = 1.0;
129  test_elem.section.y_vector() = orientation;
130  test_elem.section.init();
131  test_elem.update_residual_and_jacobian();
132 
133  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
134  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
135  }
136 
137 
138  SECTION("Internal Jacobian (stiffness matrix) is invariant to displacement solution")
139  {
140  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
141  elem_sol << 0.05727841, 0.08896581, 0.09541619, -0.03774913,
142  0.07510557, -0.07122266, -0.00979117, -0.08300009,
143  -0.03453369, -0.05487761, -0.01407677, -0.09268421;
144  test_elem.elem->set_solution(elem_sol);
145  test_elem.update_residual_and_jacobian();
146 
147  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
148  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
149  }
150 
151 
152  SECTION("Internal Jacobian (stiffness matrix) is invariant to element x-location")
153  {
154  TEST::transform_element(test_elem.mesh, coords,5.2, 0.0, 0.0,
155  1.0, 1.0, 0.0, 0.0, 0.0);
156  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
157  test_elem.update_residual_and_jacobian();
158 
159  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
160  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
161  }
162 
163 
164  SECTION("Internal Jacobian (stiffness matrix) is invariant to element y-location")
165  {
166  TEST::transform_element(test_elem.mesh, coords, 0.0, -11.5, 0.0,
167  1.0, 1.0, 0.0, 0.0, 0.0);
168  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
169  test_elem.update_residual_and_jacobian();
170 
171  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
172  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
173  }
174 
175 
176  SECTION("Internal Jacobian (stiffness matrix) is invariant to element z-location")
177  {
178  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 7.6,
179  1.0, 1.0, 0.0, 0.0, 0.0);
180  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
181  test_elem.update_residual_and_jacobian();
182 
183  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
184  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
185  }
186 
187 
188  SECTION("Internal Jacobian (stiffness matrix) checks for element aligned along y-axis")
189  {
190  /*
191  * NOTE: We could try to use the transform_element method here, but the
192  * issue is that if the sin and cos calculations are not exact, then we
193  * may not be perfectly aligned along the y axis like we want.
194  */
195  RealMatrixX new_coordinates = RealMatrixX::Zero(3,test_elem.n_nodes);
196  new_coordinates << 0.0, 0.0, -1.0, 1.0, 0.0, 0.0;
197  test_elem.update_coordinates(new_coordinates);
198  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
199  test_elem.update_residual_and_jacobian();
200 
201  // Finite difference Jacobian check
203  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
204  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
205  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
206 
207  // Symmetry check
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.transpose())));
210 
211  // Determinant check
212  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
213  }
214 
215 
216  SECTION("Internal Jacobian (stiffness matrix) checks for element aligned along z-axis")
217  {
218  /*
219  * NOTE: We could try to use the transform_element method here, but the
220  * issue is that if the sin and cos calculations are not exact, then we
221  * may not be perfectly aligned along the z axis like we want.
222  */
223  RealMatrixX new_coordinates = RealMatrixX::Zero(3,test_elem.n_nodes);
224  new_coordinates << 0.0, 0.0, 0.0, 0.0, -1.0, 1.0;
225  test_elem.update_coordinates(new_coordinates);
226  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
227  test_elem.update_residual_and_jacobian();
228 
229 
230  // Finite difference Jacobian check
232  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
233  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
234  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
235 
236  // Symmetry check
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.transpose())));
239 
240  // Determinant check
241  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
242  }
243 
244 
245  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about z-axis")
246  {
247  // Rotated 63.4 about z-axis at element's centroid
248  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
249  1.0, 1.0, 0.0, 0.0, 63.4);
250  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
251  test_elem.update_residual_and_jacobian();
252 
253  // Finite difference Jacobian check
255  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
256  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
257  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
258 
259  // Symmetry check
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.transpose())));
262 
263  // Determinant check
264  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
265  }
266 
267 
268  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about y-axis")
269  {
270  // Rotated 35.8 about y-axis at element's centroid
271  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
272  1.0, 1.0, 0.0, 35.8, 0.0);
273  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
274  test_elem.update_residual_and_jacobian();
275 
276  // Finite difference Jacobian check
278  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
279  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
280  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
281 
282  // Symmetry check
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.transpose())));
285 
286  // Determinant check
287  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
288  }
289 
290 
291  SECTION("Internal Jacobian (stiffness matrix) checks for element stretched in x-direction")
292  {
293  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
294  3.2, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0);
295  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
296  test_elem.update_residual_and_jacobian();
297 
298  // Finite difference Jacobian check
300  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
301  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
302  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
303 
304  // Symmetry check
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.transpose())));
307 
308  // Determinant check
309  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
310  }
311 
312 
313  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, and rotated")
314  {
315  // Apply arbitrary transformation to the element
316  TEST::transform_element(test_elem.mesh, coords, -5.0, 7.8, -13.1,
317  2.7, 6.4, 20.0, 47.8, -70.1);
318  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
319  test_elem.update_residual_and_jacobian();
320 
321  // Finite difference Jacobian check
323  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
324  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
325  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
326 
327  // Symmetry check
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.transpose())));
330 
331  // Determinant check
332  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
333  }
334 
335 
336  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and displaced")
337  {
338  // Apply arbitrary transformation to the element
339  TEST::transform_element(test_elem.mesh, coords, 4.1, -6.3, 7.5,
340  4.2, 1.5, -18.0, -24.8, 30.1);
341 
342  // Apply arbitrary displacement to the element
343  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
344  elem_sol << 0.08158724, 0.07991906, -0.00719128, 0.02025461,
345  -0.04602193, 0.05280159, 0.03700081, 0.04636344,
346  0.05559377, 0.06448206, 0.08919238, -0.03079122;
347  test_elem.elem->set_solution(elem_sol);
348 
349  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
350  test_elem.update_residual_and_jacobian();
351 
352  // Finite difference Jacobian check
354  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
355  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
356  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
357 
358  // Symmetry check
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.transpose())));
361 
362  // Determinant check
363  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
364  }
365 }
RealMatrixX jacobian
Matrix storage for Jacobian of the element in a perturbed/modified state.
TEST_CASE("edge2_linear_extension_bending_coupling_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
#define eps
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