MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_edge2_nonlinear_structural_thermal_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 // We need access to the protected thermal_residual method to test it
21 // NOTE: Be careful with this, it could cause unexpected problems
22 #define protected public
23 
24 // libMesh includes
25 #include "libmesh/libmesh.h"
26 #include "libmesh/elem.h"
27 #include "libmesh/dof_map.h"
28 
29 // MAST includes
30 #include "base/parameter.h"
38 #include "base/nonlinear_system.h"
39 #include "mesh/geom_elem.h"
41 
42 // Test includes
43 #include "catch.hpp"
44 #include "test_helpers.h"
46 
47 extern libMesh::LibMeshInit* p_global_init;
48 
49 
50 TEST_CASE("edge2_nonlinear_structural_thermal_jacobian",
51  "[1D],[thermoelastic],[edge],[edge2],[nonlinear],[protected]")
52 {
53  RealMatrixX coords = RealMatrixX::Zero(3, 2);
54  coords << -1.0, 1.0, 0.0,
55  0.0, 0.0, 0.0;
56  TEST::TestStructuralSingleElement1D test_elem(libMesh::EDGE2, coords);
57 
58  // Define the Uniform Temperature and Uniform Reference Temperature
59  MAST::Parameter temperature("T", 400.0);
60  MAST::Parameter ref_temperature("T0", 0.0);
61  MAST::ConstantFieldFunction temperature_f("temperature", temperature);
62  MAST::ConstantFieldFunction ref_temperature_f("ref_temperature", ref_temperature);
63 
64  // Setup the temperature change boundary condition
66  temperature_load.add(temperature_f);
67  temperature_load.add(ref_temperature_f);
68  test_elem.discipline.add_volume_load(0, temperature_load);
69 
70  // Set the strain type to non-linear for the section
72  // Set the bending operator to Euler-Bernoulli
74 
75  const Real V0 = test_elem.reference_elem->volume();
76 
77  // Calculate residual and jacobian
78  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian0, temperature_load);
79 
80  double val_margin = (test_elem.jacobian0.array().abs()).mean() * 1.490116119384766e-08;
81 
82  //libMesh::out << "R=\n" << residual << std::endl;
83  libMesh::out << "J =\n" << test_elem.jacobian0 << std::endl;
84 
85 
86  SECTION("Thermoelastic Jacobian finite difference check")
87  {
89  test_elem.elem_solution, test_elem.jacobian_fd,
90  temperature_load);
91 
92  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
93 
94  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
95  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
96  }
97 
98 
99  SECTION("Thermoelastic Jacobian should be symmetric")
100  {
101  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
102  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0.transpose())));
103  }
104 
105 
106  SECTION("Determinant of undeformed thermoelastic Jacobian should be zero")
107  {
108  REQUIRE(test_elem.jacobian0.determinant() == Approx(0.0).margin(1e-06));
109  }
110 
111 
112 // SECTION("thermal_jacobian_eigenvalue_check")
113 // {
114 // /**
115 // * Linear thermoelastic Jacobian should be independent of the
116 // * displacements and thus should be a zero matrix.
117 // */
118 // SelfAdjointEigenSolver<RealMatrixX> eigensolver(jacobian0, false);
119 // RealVectorX eigenvalues = eigensolver.eigenvalues();
120 // libMesh::out << "Eigenvalues are:\n" << eigenvalues << std::endl;
121 // uint nz = 0;
122 // for (uint i=0; i<eigenvalues.size(); i++)
123 // {
124 // if (std::abs(eigenvalues(i))<0.0001220703125)
125 // {
126 // nz++;
127 // }
128 // }
129 // REQUIRE( nz == 12);
130 // }
131 
132 
133  SECTION("Thermoelastic Jacobian is invariant to section orientation")
134  {
135  test_elem.section.clear();
136  RealVectorX orientation = RealVectorX::Zero(3);
137  orientation(2) = 1.0;
138  test_elem.section.y_vector() = orientation;
139  test_elem.section.init();
140  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
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("Thermoelastic Jacobian is invariant to displacement solution")
148  {
149  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
150  elem_sol << 0.05727841, 0.08896581, 0.09541619, -0.03774913,
151  0.07510557, -0.07122266, -0.00979117, -0.08300009,
152  -0.03453369, -0.05487761, -0.01407677, -0.09268421;
153  test_elem.elem->set_solution(elem_sol);
154  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
155 
156 
157  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
158  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
159  }
160 
161 
162  SECTION("Thermoelastic Jacobian is invariant to element x-location")
163  {
164  TEST::transform_element(test_elem.mesh, coords,5.2, 0.0, 0.0,
165  1.0, 1.0, 0.0, 0.0, 0.0);
166  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
167  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
168 
169  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
170  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
171  }
172 
173 
174  SECTION("Thermoelastic Jacobian is invariant to element y-location")
175  {
176  TEST::transform_element(test_elem.mesh, coords, 0.0, -11.5, 0.0,
177  1.0, 1.0, 0.0, 0.0, 0.0);
178  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
179  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
180 
181  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
182  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
183  }
184 
185 
186  SECTION("Thermoelastic Jacobian is invariant to element z-location")
187  {
188  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 7.6,
189  1.0, 1.0, 0.0, 0.0, 0.0);
190  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
191  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
192 
193  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
194  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
195  }
196 
197 
198  SECTION("Thermoelastic Jacobian checks for element aligned along y-axis")
199  {
200  /*
201  * NOTE: We could try to use the transform_element method here, but the
202  * issue is that if the sin and cos calculations are not exact, then we
203  * may not be perfectly aligned along the y axis like we want.
204  */
205  RealMatrixX new_coordinates = RealMatrixX::Zero(3,test_elem.n_nodes);
206  new_coordinates << 0.0, 0.0, -1.0, 1.0, 0.0, 0.0;
207  test_elem.update_coordinates(new_coordinates);
208  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
209  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
210 
211  // Finite difference Jacobian check
212  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
213  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
214  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
215  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
216 
217  // Symmetry check
218  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
219  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
220 
221  // Determinant check
222  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
223  }
224 
225  SECTION("Thermoelastic Jacobian checks for element aligned along z-axis")
226  {
227  /*
228  * NOTE: We could try to use the transform_element method here, but the
229  * issue is that if the sin and cos calculations are not exact, then we
230  * may not be perfectly aligned along the z axis like we want.
231  */
232  RealMatrixX new_coordinates = RealMatrixX::Zero(3,test_elem.n_nodes);
233  new_coordinates << 0.0, 0.0, 0.0, 0.0, -1.0, 1.0;
234  test_elem.update_coordinates(new_coordinates);
235  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
236  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
237 
238  // Finite difference Jacobian check
239  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
240  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
241  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
242  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
243 
244  // Symmetry check
245  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
246  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
247 
248  // Determinant check
249  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
250  }
251 
252 
253  SECTION("Thermoelastic Jacobian checks for element rotated about z-axis")
254  {
255  // Rotated 63.4 about z-axis at element's centroid
256  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
257  1.0, 1.0, 0.0, 0.0, 63.4);
258  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
259  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
260 
261  // Finite difference Jacobian check
262  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
263  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
264  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
265  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
266 
267  // Symmetry check
268  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
269  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
270 
271  // Determinant check
272  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
273  }
274 
275 
276  SECTION("Thermoelastic Jacobian checks for element rotated about y-axis")
277  {
278  // Rotated 35.8 about y-axis at element's centroid
279  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
280  1.0, 1.0, 0.0, 35.8, 0.0);
281  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
282  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
283 
284  // Finite difference Jacobian check
285  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
286  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
287  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
288  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
289 
290  // Symmetry check
291  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
292  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
293 
294  // Determinant check
295  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
296  }
297 
298 
299  SECTION("Thermoelastic Jacobian checks for element stretched in x-direction")
300  {
301  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
302  3.2, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0);
303  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
304  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
305 
306  // Finite difference Jacobian check
307  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
308  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
309  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
310  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
311 
312  // Symmetry check
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.transpose())));
315 
316  // Determinant check
317  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
318  }
319 
320 
321  SECTION("Thermoelastic Jacobian checks for element arbitrarily scaled, stretched, and rotated")
322  {
323  // Apply arbitrary transformation to the element
324  TEST::transform_element(test_elem.mesh, coords, -5.0, 7.8, -13.1,
325  2.7, 6.4, 20.0, 47.8, -70.1);
326  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
327  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
328 
329  // Finite difference Jacobian check
330  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
331  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
332  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
333  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
334 
335  // Symmetry check
336  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
337  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
338 
339  // Determinant check
340  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
341  }
342 
343 
344  SECTION("Thermoelastic Jacobian checks for element arbitrarily scaled, stretched, rotated, and displaced")
345  {
346  // Apply arbitrary transformation to the element
347  TEST::transform_element(test_elem.mesh, coords, 4.1, -6.3, 7.5,
348  4.2, 1.5, -18.0, -24.8, 30.1);
349 
350  // Apply arbitrary displacement to the element
351  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
352  elem_sol << 0.08158724, 0.07991906, -0.00719128, 0.02025461,
353  -0.04602193, 0.05280159, 0.03700081, 0.04636344,
354  0.05559377, 0.06448206, 0.08919238, -0.03079122;
355  test_elem.elem->set_solution(elem_sol);
356 
357  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
358  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
359 
360  // Finite difference Jacobian check
361  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
362  val_margin = (test_elem.jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
363  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
364  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
365 
366  // Symmetry check
367  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
368  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
369 
370  // Determinant check
371  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
372  }
373 }
RealMatrixX jacobian
Matrix storage for Jacobian of the element in a perturbed/modified state.
void add_volume_load(libMesh::subdomain_id_type bid, MAST::BoundaryConditionBase &load)
adds the specified volume loads for the elements with subdomain tag s_id
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
This is a scalar function whose value can be changed and one that can be used as a design variable in...
Definition: parameter.h:35
libMesh::Real Real
void set_strain(MAST::StrainType strain)
sets the type of strain to be used, which is LINEAR_STRAIN by default
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...
void add(MAST::FunctionBase &f)
adds the function to this card and returns a reference to it.
std::vector< double > eigen_matrix_to_std_vector(RealMatrixX M)
Converts an Eigen Matrix object to a std::vector.
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.
TEST_CASE("edge2_nonlinear_structural_thermal_jacobian", "[1D],[thermoelastic],[edge],[edge2],[nonlinear],[protected]")
RealVectorX residual
Vector storage for element&#39;s residual vector.
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
void approximate_thermal_jacobian_with_finite_difference(MAST::StructuralElementBase &elem, const RealVectorX &initial_elem_solution, RealMatrixX &jacobian, MAST::BoundaryConditionBase &thermal_bc)
Approximates the thermal jacobian using a 6th order accurate central finite difference scheme...
virtual bool thermal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)
Calculates the force vector and Jacobian due to thermal stresses.
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
void set_bending_model(MAST::BendingOperatorType b)
returns the bending model to be used for the 2D element
libMesh::ReplicatedMesh mesh
The actual libMesh mesh object.
Definition: mast_mesh.h:52