MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_quad4_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("quad4_linear_extension_bending_shear_structural",
46  "[quad],[quad4],[linear],[structural],[2D]")
47 {
48  RealMatrixX coords = RealMatrixX::Zero(3,4);
49  coords << -1.0, 1.0, 1.0, -1.0,
50  -1.0, -1.0, 1.0, 1.0,
51  0.0, 0.0, 0.0, 0.0;
52  TEST::TestStructuralSingleElement2D test_elem(libMesh::QUAD4, coords);
53 
54  const Real V0 = test_elem.reference_elem->volume();
55  REQUIRE(test_elem.reference_elem->volume() == TEST::get_shoelace_area(coords));
56 
57  // Set the strain type to linear for the section
59 
60  // Set the offset to zero to disable extension-bending coupling
61  test_elem.offset = 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 
69  SECTION("Internal Jacobian (stiffness matrix) finite difference check")
70  {
72  test_elem.elem_solution, test_elem.jacobian_fd);
73 
74  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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.jacobian_fd)).margin(val_margin));
78  }
79 
80 
81  SECTION("Internal Jacobian (stiffness matrix) should be symmetric")
82  {
83  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
84  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0.transpose())));
85  }
86 
87 
88  SECTION("Determinant of undeformed internal Jacobian (stiffness matrix) should be zero")
89  {
90  REQUIRE(test_elem.jacobian0.determinant() == Approx(0.0).margin(1e-06));
91  }
92 
93 
94  SECTION("Internal Jacobian (stiffness matrix) is invariant to displacement solution")
95  {
96  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
97  elem_sol << -0.04384355, 0.03969142, -0.09470648, -0.05011107,
98  -0.02989082, -0.01205296, 0.08846868, 0.04522207,
99  0.06435953, -0.07282706, 0.09307561, -0.06250143,
100  0.03332844, -0.00040089, -0.00423108, -0.07258241,
101  0.06636534, -0.08421098, -0.0705489 , -0.06004976,
102  0.03873095, -0.09194373, 0.00055061, 0.046831;
103  test_elem.elem->set_solution(elem_sol);
104  test_elem.update_residual_and_jacobian();
105 
106  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
107  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
108  }
109 
110 
111  SECTION("Internal Jacobian (stiffness matrix) is invariant to element x-location")
112  {
113  TEST::transform_element(test_elem.mesh, coords,5.2, 0.0, 0.0,
114  1.0, 1.0, 0.0, 0.0, 0.0);
115  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
116  test_elem.update_residual_and_jacobian();
117 
118  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
119  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
120  }
121 
122 
123  SECTION("Internal Jacobian (stiffness matrix) is invariant to element y-location")
124  {
125  TEST::transform_element(test_elem.mesh, coords, 0.0, -11.5, 0.0,
126  1.0, 1.0, 0.0, 0.0, 0.0);
127  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
128  test_elem.update_residual_and_jacobian();
129 
130  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
131  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
132  }
133 
134 
135  SECTION("Internal Jacobian (stiffness matrix) is invariant to element z-location")
136  {
137  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 7.6,
138  1.0, 1.0, 0.0, 0.0, 0.0);
139  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
140  test_elem.update_residual_and_jacobian();
141 
142  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
143  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
144  }
145 
146 
147  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about z-axis")
148  {
149  // Rotated 63.4 about z-axis at element's centroid
150  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
151  1.0, 1.0, 0.0, 0.0, 63.4);
152  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
153  test_elem.update_residual_and_jacobian();
154 
155  // Finite difference Jacobian check
157  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
158  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
159  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
160 
161  // Symmetry check
162  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
163  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
164 
165  // Determinant check
166  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
167  }
168 
169 
170  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about y-axis")
171  {
172  // Rotated 35.8 about y-axis at element's centroid
173  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
174  1.0, 1.0, 0.0, 35.8, 0.0);
175  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
176  test_elem.update_residual_and_jacobian();
177 
178  // Finite difference Jacobian check
180  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
181  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
182  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
183 
184  // Symmetry check
185  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
186  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
187 
188  // Determinant check
189  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
190  }
191 
192 
193  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about x-axis")
194  {
195  // Rotated 15.8 about x-axis at element's centroid
196  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
197  1.0, 1.0, 15.8, 0.0, 0.0);
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 sheared in x-direction\"")
217  {
218  // Shear element in x-direction.
219  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
220  1.0, 1.0, 0.0, 0.0, 0.0, 6.7, 0.0);
221  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
222  test_elem.update_residual_and_jacobian();
223 
224  // Finite difference Jacobian check
226  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
227  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
228  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
229 
230  // Symmetry check
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.transpose())));
233 
234  // Determinant check
235  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
236  }
237 
238 
239  SECTION("\"Internal Jacobian (stiffness matrix) checks for element sheared in y-direction\"")
240  {
241  // Shear element in x-direction.
242  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
243  1.0, 1.0, 0.0, 0.0, 0.0, 0.0, -11.2);
244  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
245  test_elem.update_residual_and_jacobian();
246 
247  // Finite difference Jacobian check
249  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
250  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
251  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
252 
253  // Symmetry check
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.transpose())));
256 
257  // Determinant check
258  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
259  }
260 
261 
262  SECTION("Internal Jacobian (stiffness matrix) checks for element stretched in x-direction")
263  {
264  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
265  3.2, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0);
266  REQUIRE_FALSE(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  SECTION("Internal Jacobian (stiffness matrix) checks for element stretched in y-direction")
284  {
285  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
286  1.0, 0.64, 0.0, 0.0, 0.0, 0.0, 0.0);
287  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
288  test_elem.update_residual_and_jacobian();
289 
290  // Finite difference Jacobian check
292  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
293  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
294  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
295 
296  // Symmetry check
297  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
298  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
299 
300  // Determinant check
301  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
302  }
303 
304 
305  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and sheared")
306  {
307  // Apply arbitrary transformation to the element
308  TEST::transform_element(test_elem.mesh, coords, -5.0, 7.8, -13.1,
309  2.7, 6.4, 20.0, 47.8, -70.1, 5.7, -6.3);
310  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
311  test_elem.update_residual_and_jacobian();
312 
313  // Finite difference Jacobian check
315  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
316  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
317  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
318 
319  // Symmetry check
320  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
321  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
322 
323  // Determinant check
324  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
325  }
326 
327 
328  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and displaced")
329  {
330  // Apply arbitrary transformation to the element
331  TEST::transform_element(test_elem.mesh, coords, 4.1, -6.3, 7.5, 4.2, 1.5, -18.0, -24.8,
332  30.1, -3.2, 5.4);
333 
334  // Apply arbitrary displacement to the element
335  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
336  elem_sol << -0.04384355, 0.03969142, -0.09470648, -0.05011107,
337  -0.02989082, -0.01205296, 0.08846868, 0.04522207,
338  0.06435953, -0.07282706, 0.09307561, -0.06250143,
339  0.03332844, -0.00040089, -0.00423108, -0.07258241,
340  0.06636534, -0.08421098, -0.0705489 , -0.06004976,
341  0.03873095, -0.09194373, 0.00055061, 0.046831;
342  test_elem.elem->set_solution(elem_sol);
343 
344  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
345  test_elem.update_residual_and_jacobian();
346 
347  // Finite difference Jacobian check
349  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
352 
353  // Symmetry check
354  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
355  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
356 
357  // Determinant check
358  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
359  }
360 }
MAST::Parameter offset
Section offset.
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
Real get_shoelace_area(RealMatrixX X)
Calcualtes the area of a 2D polygon using the shoelace formula.
MAST::Solid2DSectionElementPropertyCard section
libMesh::Real Real
void set_strain(MAST::StrainType strain)
sets the type of strain to be used, which is LINEAR_STRAIN by default
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...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
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...
RealMatrixX jacobian0
Matrix storage for Jacobian of baseline/undeformed element.
Matrix< Real, Dynamic, 1 > RealVectorX
libMesh::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26
TEST_CASE("quad4_linear_extension_bending_shear_structural", "[quad],[quad4],[linear],[structural],[2D]")
libMesh::ReplicatedMesh mesh
The actual libMesh mesh object.
Definition: mast_mesh.h:52
RealMatrixX jacobian_fd
Matrix storage for element Jacobian approximated by finite difference.