MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_quad4_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("quad4_linear_extension_bending_structural",
46  "[quad],[quad4],[linear],[structural],[2D],[element]")
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 shear coefficient to zero to disable transverse shear stiffness
61  test_elem.kappa = 0.0;
62 
63  // Set the offset to zero to disable extension-bending coupling
64  test_elem.offset = 0.0;
65 
66  // Update residual and Jacobian storage since we have modified baseline test element properties.
68 
69  double val_margin = (test_elem.jacobian0.array().abs()).mean() * 1.490116119384766e-08;
70 
71 
72  SECTION("Internal Jacobian (stiffness matrix) finite difference check")
73  {
75  test_elem.elem_solution, test_elem.jacobian_fd);
76 
77  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
78 
79  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
80  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
81  }
82 
83 
84  SECTION("Internal Jacobian (stiffness matrix) should be symmetric")
85  {
86  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
87  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0.transpose())));
88  }
89 
90 
91  SECTION("Determinant of undeformed internal Jacobian (stiffness matrix) should be zero")
92  {
93  REQUIRE(test_elem.jacobian0.determinant() == Approx(0.0).margin(1e-06));
94  }
95 
96 
97  SECTION("Internal Jacobian (stiffness matrix) is invariant to displacement solution")
98  {
99  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
100  elem_sol << -0.04384355, 0.03969142, -0.09470648, -0.05011107,
101  -0.02989082, -0.01205296, 0.08846868, 0.04522207,
102  0.06435953, -0.07282706, 0.09307561, -0.06250143,
103  0.03332844, -0.00040089, -0.00423108, -0.07258241,
104  0.06636534, -0.08421098, -0.0705489 , -0.06004976,
105  0.03873095, -0.09194373, 0.00055061, 0.046831;
106  test_elem.elem->set_solution(elem_sol);
107  test_elem.update_residual_and_jacobian();
108 
109  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
110  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
111  }
112 
113 
114  SECTION("Internal Jacobian (stiffness matrix) is invariant to element x-location")
115  {
116  TEST::transform_element(test_elem.mesh, coords,5.2, 0.0, 0.0,
117  1.0, 1.0, 0.0, 0.0, 0.0);
118  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
119  test_elem.update_residual_and_jacobian();
120 
121  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
122  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
123  }
124 
125 
126  SECTION("Internal Jacobian (stiffness matrix) is invariant to element y-location")
127  {
128  TEST::transform_element(test_elem.mesh, coords, 0.0, -11.5, 0.0,
129  1.0, 1.0, 0.0, 0.0, 0.0);
130  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
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 element z-location")
139  {
140  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 7.6,
141  1.0, 1.0, 0.0, 0.0, 0.0);
142  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
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) checks for element rotated about z-axis")
151  {
152  // Rotated 63.4 about z-axis at element's centroid
153  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
154  1.0, 1.0, 0.0, 0.0, 63.4);
155  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
156  test_elem.update_residual_and_jacobian();
157 
158  // Finite difference Jacobian check
160  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
161  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
162  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
163 
164  // Symmetry check
165  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
166  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
167 
168  // Determinant check
169  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
170  }
171 
172 
173  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about y-axis")
174  {
175  // Rotated 35.8 about y-axis at element's centroid
176  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
177  1.0, 1.0, 0.0, 35.8, 0.0);
178  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
179  test_elem.update_residual_and_jacobian();
180 
181  // Finite difference Jacobian check
183  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
184  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
185  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
186 
187  // Symmetry check
188  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
189  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
190 
191  // Determinant check
192  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
193  }
194 
195 
196  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about x-axis")
197  {
198  // Rotated 15.8 about x-axis at element's centroid
199  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
200  1.0, 1.0, 15.8, 0.0, 0.0);
201  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
202  test_elem.update_residual_and_jacobian();
203 
204  // Finite difference Jacobian check
206  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
207  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
208  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
209 
210  // Symmetry check
211  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
212  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
213 
214  // Determinant check
215  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
216  }
217 
218 
219  SECTION("\"Internal Jacobian (stiffness matrix) checks for element sheared in x-direction\"")
220  {
221  // Shear element in x-direction.
222  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
223  1.0, 1.0, 0.0, 0.0, 0.0, 6.7, 0.0);
224  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
225  test_elem.update_residual_and_jacobian();
226 
227  // Finite difference Jacobian check
229  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
230  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
231  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
232 
233  // Symmetry check
234  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
235  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
236 
237  // Determinant check
238  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
239  }
240 
241 
242  SECTION("\"Internal Jacobian (stiffness matrix) checks for element sheared in y-direction\"")
243  {
244  // Shear element in x-direction.
245  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
246  1.0, 1.0, 0.0, 0.0, 0.0, 0.0, -11.2);
247  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
248  test_elem.update_residual_and_jacobian();
249 
250  // Finite difference Jacobian check
252  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
253  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
254  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
255 
256  // Symmetry check
257  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
258  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
259 
260  // Determinant check
261  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
262  }
263 
264 
265  SECTION("Internal Jacobian (stiffness matrix) checks for element stretched in x-direction")
266  {
267  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
268  3.2, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0);
269  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
270  test_elem.update_residual_and_jacobian();
271 
272  // Finite difference Jacobian check
274  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
275  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
276  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
277 
278  // Symmetry check
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.transpose())));
281 
282  // Determinant check
283  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
284  }
285 
286  SECTION("Internal Jacobian (stiffness matrix) checks for element stretched in y-direction")
287  {
288  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
289  1.0, 0.64, 0.0, 0.0, 0.0, 0.0, 0.0);
290  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
291  test_elem.update_residual_and_jacobian();
292 
293  // Finite difference Jacobian check
295  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
298 
299  // Symmetry check
300  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
301  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
302 
303  // Determinant check
304  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
305  }
306 
307 
308  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and sheared")
309  {
310  // Apply arbitrary transformation to the element
311  TEST::transform_element(test_elem.mesh, coords, -5.0, 7.8, -13.1,
312  2.7, 6.4, 20.0, 47.8, -70.1, 5.7, -6.3);
313  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
314  test_elem.update_residual_and_jacobian();
315 
316  // Finite difference Jacobian check
318  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
321 
322  // Symmetry check
323  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
324  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
325 
326  // Determinant check
327  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
328  }
329 
330 
331  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and displaced")
332  {
333  // Apply arbitrary transformation to the element
334  TEST::transform_element(test_elem.mesh, coords, 4.1, -6.3, 7.5, 4.2, 1.5, -18.0, -24.8,
335  30.1, -3.2, 5.4);
336 
337  // Apply arbitrary displacement to the element
338  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
339  elem_sol << -0.04384355, 0.03969142, -0.09470648, -0.05011107,
340  -0.02989082, -0.01205296, 0.08846868, 0.04522207,
341  0.06435953, -0.07282706, 0.09307561, -0.06250143,
342  0.03332844, -0.00040089, -0.00423108, -0.07258241,
343  0.06636534, -0.08421098, -0.0705489 , -0.06004976,
344  0.03873095, -0.09194373, 0.00055061, 0.046831;
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 }
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
libMesh::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26
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
MAST::Parameter kappa
Shear coefficient.
TEST_CASE("quad4_linear_extension_bending_structural", "[quad],[quad4],[linear],[structural],[2D],[element]")
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.