MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_quad4_linear_structural_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_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  // Update residual and Jacobian storage since we have modified baseline test element properties.
62 
63  double val_margin = (test_elem.jacobian0.array().abs()).mean() * 1.490116119384766e-08;
64 
65 
66  SECTION("Internal Jacobian (stiffness matrix) finite difference check")
67  {
69  test_elem.elem_solution, test_elem.jacobian_fd);
70 
71  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
72 
73  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
74  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
75  }
76 
77 
78  SECTION("Internal Jacobian (stiffness matrix) should be symmetric")
79  {
80  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
81  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0.transpose())));
82  }
83 
84 
85  SECTION("Determinant of undeformed internal Jacobian (stiffness matrix) should be zero")
86  {
87  REQUIRE(test_elem.jacobian0.determinant() == Approx(0.0).margin(1e-06));
88  }
89 
90 
91  SECTION("Internal Jacobian (stiffness matrix) is invariant to displacement solution")
92  {
93  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
94  elem_sol << -0.04384355, 0.03969142, -0.09470648, -0.05011107,
95  -0.02989082, -0.01205296, 0.08846868, 0.04522207,
96  0.06435953, -0.07282706, 0.09307561, -0.06250143,
97  0.03332844, -0.00040089, -0.00423108, -0.07258241,
98  0.06636534, -0.08421098, -0.0705489 , -0.06004976,
99  0.03873095, -0.09194373, 0.00055061, 0.046831;
100  test_elem.elem->set_solution(elem_sol);
101  test_elem.update_residual_and_jacobian();
102 
103  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
104  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
105  }
106 
107 
108  SECTION("Internal Jacobian (stiffness matrix) is invariant to element x-location")
109  {
110  TEST::transform_element(test_elem.mesh, coords,5.2, 0.0, 0.0,
111  1.0, 1.0, 0.0, 0.0, 0.0);
112  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
113  test_elem.update_residual_and_jacobian();
114 
115  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
116  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
117  }
118 
119 
120  SECTION("Internal Jacobian (stiffness matrix) is invariant to element y-location")
121  {
122  TEST::transform_element(test_elem.mesh, coords, 0.0, -11.5, 0.0,
123  1.0, 1.0, 0.0, 0.0, 0.0);
124  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
125  test_elem.update_residual_and_jacobian();
126 
127  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
128  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
129  }
130 
131 
132  SECTION("Internal Jacobian (stiffness matrix) is invariant to element z-location")
133  {
134  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 7.6,
135  1.0, 1.0, 0.0, 0.0, 0.0);
136  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
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) checks for element rotated about z-axis")
145  {
146  // Rotated 63.4 about z-axis at element's centroid
147  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
148  1.0, 1.0, 0.0, 0.0, 63.4);
149  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
150  test_elem.update_residual_and_jacobian();
151 
152  // Finite difference Jacobian check
154  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
155  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
156  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
157 
158  // Symmetry check
159  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
160  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
161 
162  // Determinant check
163  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
164  }
165 
166 
167  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about y-axis")
168  {
169  // Rotated 35.8 about y-axis at element's centroid
170  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
171  1.0, 1.0, 0.0, 35.8, 0.0);
172  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
173  test_elem.update_residual_and_jacobian();
174 
175  // Finite difference Jacobian check
177  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
178  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
179  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
180 
181  // Symmetry check
182  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
183  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
184 
185  // Determinant check
186  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
187  }
188 
189 
190  SECTION("Internal Jacobian (stiffness matrix) checks for element rotated about x-axis")
191  {
192  // Rotated 15.8 about x-axis at element's centroid
193  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
194  1.0, 1.0, 15.8, 0.0, 0.0);
195  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
196  test_elem.update_residual_and_jacobian();
197 
198  // Finite difference Jacobian check
200  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
201  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
202  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
203 
204  // Symmetry check
205  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
206  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
207 
208  // Determinant check
209  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
210  }
211 
212 
213  SECTION("\"Internal Jacobian (stiffness matrix) checks for element sheared in x-direction\"")
214  {
215  // Shear element in x-direction.
216  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
217  1.0, 1.0, 0.0, 0.0, 0.0, 6.7, 0.0);
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 sheared in y-direction\"")
237  {
238  // Shear element in x-direction.
239  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
240  1.0, 1.0, 0.0, 0.0, 0.0, 0.0, -11.2);
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 stretched in x-direction")
260  {
261  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
262  3.2, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0);
263  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
264  test_elem.update_residual_and_jacobian();
265 
266  // Finite difference Jacobian check
268  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
269  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
270  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
271 
272  // Symmetry check
273  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
274  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
275 
276  // Determinant check
277  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
278  }
279 
280  SECTION("Internal Jacobian (stiffness matrix) checks for element stretched in y-direction")
281  {
282  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
283  1.0, 0.64, 0.0, 0.0, 0.0, 0.0, 0.0);
284  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
285  test_elem.update_residual_and_jacobian();
286 
287  // Finite difference Jacobian check
289  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
290  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
291  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
292 
293  // Symmetry check
294  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
295  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
296 
297  // Determinant check
298  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
299  }
300 
301 
302  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and sheared")
303  {
304  // Apply arbitrary transformation to the element
305  TEST::transform_element(test_elem.mesh, coords, -5.0, 7.8, -13.1,
306  2.7, 6.4, 20.0, 47.8, -70.1, 5.7, -6.3);
307  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
308  test_elem.update_residual_and_jacobian();
309 
310  // Finite difference Jacobian check
312  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
313  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
314  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
315 
316  // Symmetry check
317  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
318  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
319 
320  // Determinant check
321  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
322  }
323 
324 
325  SECTION("Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and displaced")
326  {
327  // Apply arbitrary transformation to the element
328  TEST::transform_element(test_elem.mesh, coords, 4.1, -6.3, 7.5, 4.2, 1.5, -18.0, -24.8,
329  30.1, -3.2, 5.4);
330 
331  // Apply arbitrary displacement to the element
332  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
333  elem_sol << -0.04384355, 0.03969142, -0.09470648, -0.05011107,
334  -0.02989082, -0.01205296, 0.08846868, 0.04522207,
335  0.06435953, -0.07282706, 0.09307561, -0.06250143,
336  0.03332844, -0.00040089, -0.00423108, -0.07258241,
337  0.06636534, -0.08421098, -0.0705489 , -0.06004976,
338  0.03873095, -0.09194373, 0.00055061, 0.046831;
339  test_elem.elem->set_solution(elem_sol);
340 
341  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
342  test_elem.update_residual_and_jacobian();
343 
344  // Finite difference Jacobian check
346  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
347  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
348  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
349 
350  // Symmetry check
351  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
352  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
353 
354  // Determinant check
355  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
356  }
357 }
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("quad4_linear_structural", "[quad],[quad4],[linear],[structural],[2D],[element]")
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...
libMesh::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26
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::ReplicatedMesh mesh
The actual libMesh mesh object.
Definition: mast_mesh.h:52
RealMatrixX jacobian_fd
Matrix storage for element Jacobian approximated by finite difference.