MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_edge2_linear_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_linear_structural_thermal_jacobian",
51  "[1D],[thermoelastic],[edge],[edge2],[linear],[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  const Real V0 = test_elem.reference_elem->volume();
71 
72  // Calculate residual and jacobian
73  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian0, temperature_load);
74 
75  double val_margin = (test_elem.jacobian0.array().abs()).mean() * 1.490116119384766e-08;
76 
77  libMesh::out << "R=\n" << test_elem.residual << std::endl;
78  libMesh::out << "J =\n" << test_elem.jacobian0 << std::endl;
79 
80  SECTION("thermal_jacobian_is_zero_matrix")
81  {
82  // Approximate Jacobian with Finite Difference
83  RealMatrixX zero_matrix = RealMatrixX::Zero(test_elem.n_dofs, test_elem.n_dofs);
84 
85  std::vector<double> test = TEST::eigen_matrix_to_std_vector(test_elem.jacobian0);
86  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(zero_matrix);
87 
88  REQUIRE_THAT( test, Catch::Approx<double>(truth).margin(0.0) );
89  }
90 
91 
92  SECTION("Thermoelastic Jacobian finite difference check")
93  {
95  test_elem.elem_solution, test_elem.jacobian_fd,
96  temperature_load);
97 
98  val_margin = 1.0e-5;
99 
100  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
101  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
102  }
103 
104 
105  SECTION("Thermoelastic Jacobian should be symmetric")
106  {
107  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0),
108  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0.transpose())));
109  }
110 
111 
112  SECTION("Determinant of undeformed thermoelastic Jacobian should be zero")
113  {
114  REQUIRE(test_elem.jacobian0.determinant() == Approx(0.0).margin(1e-06));
115  }
116 
117 
118  SECTION("Thermoelastic Jacobian eigenvalue check")
119  {
120  /*
121  * Linear thermoelastic Jacobian should be independent of the
122  * displacements and thus should be a zero matrix.
123  */
124  SelfAdjointEigenSolver<RealMatrixX> eigensolver(test_elem.jacobian0, false);
125  RealVectorX eigenvalues = eigensolver.eigenvalues();
126  libMesh::out << "Eigenvalues are:\n" << eigenvalues << std::endl;
127  uint nz = 0;
128  for (uint i=0; i<eigenvalues.size(); i++)
129  {
130  if (std::abs(eigenvalues(i))<0.0001220703125)
131  {
132  nz++;
133  }
134  }
135  REQUIRE( nz == 12);
136  }
137 
138 
139  SECTION("Thermoelastic Jacobian is invariant to section orientation")
140  {
141  test_elem.section.clear();
142  RealVectorX orientation = RealVectorX::Zero(3);
143  orientation(2) = 1.0;
144  test_elem.section.y_vector() = orientation;
145  test_elem.section.init();
146  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
147 
148  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
149  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
150  }
151 
152 
153  SECTION("Thermoelastic Jacobian is invariant to displacement solution")
154  {
155  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
156  elem_sol << 0.05727841, 0.08896581, 0.09541619, -0.03774913,
157  0.07510557, -0.07122266, -0.00979117, -0.08300009,
158  -0.03453369, -0.05487761, -0.01407677, -0.09268421;
159  test_elem.elem->set_solution(elem_sol);
160  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
161 
162 
163  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
164  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
165  }
166 
167 
168  SECTION("Thermoelastic Jacobian is invariant to element x-location")
169  {
170  TEST::transform_element(test_elem.mesh, coords,5.2, 0.0, 0.0,
171  1.0, 1.0, 0.0, 0.0, 0.0);
172  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
173  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
174 
175  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
176  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
177  }
178 
179 
180  SECTION("Thermoelastic Jacobian is invariant to element y-location")
181  {
182  TEST::transform_element(test_elem.mesh, coords, 0.0, -11.5, 0.0,
183  1.0, 1.0, 0.0, 0.0, 0.0);
184  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
185  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
186 
187  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
188  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
189  }
190 
191 
192  SECTION("Thermoelastic Jacobian is invariant to element z-location")
193  {
194  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 7.6,
195  1.0, 1.0, 0.0, 0.0, 0.0);
196  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
197  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
198 
199  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
200  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian0)).margin(val_margin));
201  }
202 
203 
204  SECTION("Thermoelastic Jacobian checks for element aligned along y-axis")
205  {
206  /*
207  * NOTE: We could try to use the transform_element method here, but the
208  * issue is that if the sin and cos calculations are not exact, then we
209  * may not be perfectly aligned along the y axis like we want.
210  */
211  RealMatrixX new_coordinates = RealMatrixX::Zero(3,test_elem.n_nodes);
212  new_coordinates << 0.0, 0.0, -1.0, 1.0, 0.0, 0.0;
213  test_elem.update_coordinates(new_coordinates);
214  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
215  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
216 
217  // Finite difference Jacobian check
218  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
219  val_margin = 1.0e-5;
220  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
221  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
222 
223  // Symmetry check
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.transpose())));
226 
227  // Determinant check
228  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
229  }
230 
231  SECTION("Thermoelastic Jacobian checks for element aligned along z-axis")
232  {
233  /*
234  * NOTE: We could try to use the transform_element method here, but the
235  * issue is that if the sin and cos calculations are not exact, then we
236  * may not be perfectly aligned along the z axis like we want.
237  */
238  RealMatrixX new_coordinates = RealMatrixX::Zero(3,test_elem.n_nodes);
239  new_coordinates << 0.0, 0.0, 0.0, 0.0, -1.0, 1.0;
240  test_elem.update_coordinates(new_coordinates);
241  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
242  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
243 
244  // Finite difference Jacobian check
245  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
246  val_margin = 1.0e-5;
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  SECTION("Thermoelastic Jacobian checks for element rotated about z-axis")
258  {
259  // Rotated 63.4 about z-axis at element's centroid
260  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
261  1.0, 1.0, 0.0, 0.0, 63.4);
262  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
263  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
264 
265  // Finite difference Jacobian check
266  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
267  val_margin = 1.0e-6;
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_fd)).margin(val_margin));
270 
271  // Symmetry check
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.transpose())));
274 
275  // Determinant check
276  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
277  }
278 
279 
280  SECTION("Thermoelastic Jacobian checks for element rotated about y-axis")
281  {
282  // Rotated 35.8 about y-axis at element's centroid
283  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
284  1.0, 1.0, 0.0, 35.8, 0.0);
285  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
286  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
287 
288  // Finite difference Jacobian check
289  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
290  val_margin = 1.0e-6;
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_fd)).margin(val_margin));
293 
294  // Symmetry check
295  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
296  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
297 
298  // Determinant check
299  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
300  }
301  SECTION("Thermoelastic Jacobian checks for element stretched in x-direction")
302  {
303  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0,
304  3.2, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0);
305  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
306  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
307 
308  // Finite difference Jacobian check
309  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
310  val_margin = 1.0e-6;
311  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
312  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
313 
314  // Symmetry check
315  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
316  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
317 
318  // Determinant check
319  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
320  }
321 
322 
323  SECTION("Thermoelastic Jacobian checks for element arbitrarily scaled, stretched, and rotated")
324  {
325  // Apply arbitrary transformation to the element
326  TEST::transform_element(test_elem.mesh, coords, -5.0, 7.8, -13.1,
327  2.7, 6.4, 20.0, 47.8, -70.1);
328  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
329  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
330 
331  // Finite difference Jacobian check
332  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
333  val_margin = 1.0e-6;
334  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
335  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
336 
337  // Symmetry check
338  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
339  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
340 
341  // Determinant check
342  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
343  }
344 
345 
346  SECTION("Thermoelastic Jacobian checks for element arbitrarily scaled, stretched, rotated, and displaced")
347  {
348  // Apply arbitrary transformation to the element
349  TEST::transform_element(test_elem.mesh, coords, 4.1, -6.3, 7.5,
350  4.2, 1.5, -18.0, -24.8, 30.1);
351 
352  // Apply arbitrary displacement to the element
353  RealVectorX elem_sol = RealVectorX::Zero(test_elem.n_dofs);
354  elem_sol << 0.08158724, 0.07991906, -0.00719128, 0.02025461,
355  -0.04602193, 0.05280159, 0.03700081, 0.04636344,
356  0.05559377, 0.06448206, 0.08919238, -0.03079122;
357  test_elem.elem->set_solution(elem_sol);
358 
359  REQUIRE_FALSE(test_elem.reference_elem->volume() == Approx(V0));
360  test_elem.elem->thermal_residual(true, test_elem.residual, test_elem.jacobian, temperature_load);
361 
362  // Finite difference Jacobian check
363  TEST::approximate_thermal_jacobian_with_finite_difference(*test_elem.elem, test_elem.elem_solution, test_elem.jacobian_fd, temperature_load);
364  val_margin = 1.0e-6;
365  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
366  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian_fd)).margin(val_margin));
367 
368  // Symmetry check
369  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(test_elem.jacobian),
370  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.jacobian.transpose())));
371 
372  // Determinant check
373  REQUIRE(test_elem.jacobian.determinant() == Approx(0.0).margin(1e-06));
374  }
375 }
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::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26
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
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_linear_structural_thermal_jacobian", "[1D],[thermoelastic],[edge],[edge2],[linear],[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.
void update_coordinates(RealMatrixX &new_coordinates)
Update the nodal coordinates in the mesh.
Definition: mast_mesh.h:113
libMesh::ReplicatedMesh mesh
The actual libMesh mesh object.
Definition: mast_mesh.h:52