MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_quad4_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("quad4_linear_extension_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 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  // Set the bending operator to no_bending to disable bending stiffness. This also disables
67  // transverse shear stiffness.
69 
70  // Update residual and Jacobian storage since we have modified baseline test element properties.
72 
73  double val_margin = (test_elem.jacobian0.array().abs()).mean() * 1.490116119384766e-08;
74 
75 
76  SECTION("Internal Jacobian (stiffness matrix) finite difference check")
77  {
79  test_elem.elem_solution, test_elem.jacobian_fd);
80 
81  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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.jacobian_fd)).margin(val_margin));
85  }
86 
87 
88  SECTION("Internal Jacobian (stiffness matrix) should be symmetric")
89  {
90  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
91  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0.transpose())));
92  }
93 
94 
95  SECTION("Determinant of undeformed internal Jacobian (stiffness matrix) should be zero")
96  {
97  REQUIRE(test_elem.jacobian0.determinant() == Approx(0.0).margin(1e-06));
98  }
99 
100 
101  SECTION("Internal Jacobian (stiffness matrix) is invariant to displacement solution")
102  {
103  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
104  elem_sol << -0.04384355, 0.03969142, -0.09470648, -0.05011107,
105  -0.02989082, -0.01205296, 0.08846868, 0.04522207,
106  0.06435953, -0.07282706, 0.09307561, -0.06250143,
107  0.03332844, -0.00040089, -0.00423108, -0.07258241,
108  0.06636534, -0.08421098, -0.0705489 , -0.06004976,
109  0.03873095, -0.09194373, 0.00055061, 0.046831;
110  test_elem.elem->set_solution(elem_sol);
111  test_elem.update_residual_and_jacobian();
112 
113  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
114  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
115  }
116 
117 
118  SECTION("Internal Jacobian (stiffness matrix) is invariant to element x-location")
119  {
120  TEST::transform_element(test_elem.mesh, coords,5.2, 0.0, 0.0,
121  1.0, 1.0, 0.0, 0.0, 0.0);
122  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
123  test_elem.update_residual_and_jacobian();
124 
125  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
126  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
127  }
128 
129 
130  SECTION("Internal Jacobian (stiffness matrix) is invariant to element y-location")
131  {
132  TEST::transform_element(test_elem.mesh, coords, 0.0, -11.5, 0.0,
133  1.0, 1.0, 0.0, 0.0, 0.0);
134  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
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 element z-location")
143  {
144  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 7.6,
145  1.0, 1.0, 0.0, 0.0, 0.0);
146  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
147  test_elem.update_residual_and_jacobian();
148 
149  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
150  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
151  }
152 
153 
154  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about z-axis")
155  {
156  // Rotated 63.4 about z-axis at element's centroid
157  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
158  1.0, 1.0, 0.0, 0.0, 63.4);
159  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
160  test_elem.update_residual_and_jacobian();
161 
162  // Finite difference Jacobian check
164  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
167 
168  // Symmetry check
169  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
170  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
171 
172  // Determinant check
173  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
174  }
175 
176 
177  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about y-axis")
178  {
179  // Rotated 35.8 about y-axis at element's centroid
180  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
181  1.0, 1.0, 0.0, 35.8, 0.0);
182  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
183  test_elem.update_residual_and_jacobian();
184 
185  // Finite difference Jacobian check
187  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
190 
191  // Symmetry check
192  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
193  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
194 
195  // Determinant check
196  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
197  }
198 
199 
200  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about x-axis")
201  {
202  // Rotated 15.8 about x-axis at element's centroid
203  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
204  1.0, 1.0, 15.8, 0.0, 0.0);
205  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
206  test_elem.update_residual_and_jacobian();
207 
208  // Finite difference Jacobian check
210  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
213 
214  // Symmetry check
215  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
216  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
217 
218  // Determinant check
219  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
220  }
221 
222 
223  SECTION("\"Internal Jacobian (stiffness matrix) checks for element sheared in x-direction\"")
224  {
225  // Shear element in x-direction.
226  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
227  1.0, 1.0, 0.0, 0.0, 0.0, 6.7, 0.0);
228  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
229  test_elem.update_residual_and_jacobian();
230 
231  // Finite difference Jacobian check
233  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
236 
237  // Symmetry check
238  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
239  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
240 
241  // Determinant check
242  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
243  }
244 
245 
246  SECTION("\"Internal Jacobian (stiffness matrix) checks for element sheared in y-direction\"")
247  {
248  // Shear element in x-direction.
249  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
250  1.0, 1.0, 0.0, 0.0, 0.0, 0.0, -11.2);
251  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
252  test_elem.update_residual_and_jacobian();
253 
254  // Finite difference Jacobian check
256  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
259 
260  // Symmetry check
261  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
262  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
263 
264  // Determinant check
265  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
266  }
267 
268 
269  SECTION("Internal Jacobian (stiffness matrix) checks for element stretched in x-direction")
270  {
271  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
272  3.2, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0);
273  REQUIRE_FALSE(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  SECTION("Internal Jacobian (stiffness matrix) checks for element stretched in y-direction")
291  {
292  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
293  1.0, 0.64, 0.0, 0.0, 0.0, 0.0, 0.0);
294  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
295  test_elem.update_residual_and_jacobian();
296 
297  // Finite difference Jacobian check
299  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
302 
303  // Symmetry check
304  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
305  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
306 
307  // Determinant check
308  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
309  }
310 
311 
312  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and sheared")
313  {
314  // Apply arbitrary transformation to the element
315  TEST::transform_element(test_elem.mesh, coords, -5.0, 7.8, -13.1,
316  2.7, 6.4, 20.0, 47.8, -70.1, 5.7, -6.3);
317  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
318  test_elem.update_residual_and_jacobian();
319 
320  // Finite difference Jacobian check
322  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
325 
326  // Symmetry check
327  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
328  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
329 
330  // Determinant check
331  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
332  }
333 
334 
335  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and displaced")
336  {
337  // Apply arbitrary transformation to the element
338  TEST::transform_element(test_elem.mesh, coords, 4.1, -6.3, 7.5, 4.2, 1.5, -18.0, -24.8,
339  30.1, -3.2, 5.4);
340 
341  // Apply arbitrary displacement to the element
342  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
343  elem_sol << -0.04384355, 0.03969142, -0.09470648, -0.05011107,
344  -0.02989082, -0.01205296, 0.08846868, 0.04522207,
345  0.06435953, -0.07282706, 0.09307561, -0.06250143,
346  0.03332844, -0.00040089, -0.00423108, -0.07258241,
347  0.06636534, -0.08421098, -0.0705489 , -0.06004976,
348  0.03873095, -0.09194373, 0.00055061, 0.046831;
349  test_elem.elem->set_solution(elem_sol);
350 
351  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
352  test_elem.update_residual_and_jacobian();
353 
354  // Finite difference Jacobian check
356  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
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_fd)).margin(val_margin));
359 
360  // Symmetry check
361  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
362  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
363 
364  // Determinant check
365  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
366  }
367 }
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...
TEST_CASE("quad4_linear_extension_structural", "[quad],[quad4],[linear],[structural],[2D]")
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
void set_bending_model(MAST::BendingOperatorType b)
returns the bending model to be used for the 2D element
MAST::Parameter kappa
Shear coefficient.
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.