MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_structural_element_1d.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 TEST_CASE("structural_element_1d_base_tests",
45  "[1D],[structural],[base]")
46 {
47  RealMatrixX coords = RealMatrixX::Zero(3, 2);
48  coords << -1.0, 1.0, 0.0,
49  0.0, 0.0, 0.0;
50  TEST::TestStructuralSingleElement1D test_elem(libMesh::EDGE2, coords);
51 
52  SECTION("Element returns proper number of strain components")
53  {
54  REQUIRE(test_elem.elem->n_direct_strain_components() == 2);
55  REQUIRE(test_elem.elem->n_von_karman_strain_components() == 2);
56  }
57 
58  SECTION("Incompatible modes flag returns false for basic element")
59  {
60  REQUIRE_FALSE(test_elem.elem->if_incompatible_modes() );
61  }
62 
63  SECTION("Isotropic flag returns true for basic element")
64  {
65  const MAST::ElementPropertyCardBase& elem_section = test_elem.elem->elem_property();
66  CHECK( elem_section.if_isotropic() );
67  }
68 
69  SECTION("Check setting/getting element local solution")
70  {
71  const libMesh::DofMap& dof_map = test_elem.assembly.system().get_dof_map();
72  std::vector<libMesh::dof_id_type> dof_indices;
73  dof_map.dof_indices (test_elem.reference_elem, dof_indices);
74  uint n_dofs = uint(dof_indices.size());
75 
76  RealVectorX elem_solution = 5.3*RealVectorX::Ones(n_dofs);
77  test_elem.elem->set_solution(elem_solution);
78 
79  const RealVectorX& local_solution = test_elem.elem->local_solution();
80 
81  // Convert the test and truth Eigen::Matrix objects to std::vector
82  // since Catch2 has built in methods to compare vectors
83  std::vector<double> test = TEST::eigen_matrix_to_std_vector(elem_solution);
84  std::vector<double> truth = TEST::eigen_matrix_to_std_vector(local_solution);
85 
86  // Floating point approximations are diffcult to compare since the
87  // values typically aren't exactly equal due to numerical error.
88  // Therefore, we use the Approx comparison instead of Equals
89  REQUIRE_THAT( test, Catch::Approx<double>(truth) );
90  }
91 
92  SECTION("Check setting/getting element local sensitivity solution")
93  {
94  const libMesh::DofMap& dof_map = test_elem.assembly.system().get_dof_map();
95  std::vector<libMesh::dof_id_type> dof_indices;
96  dof_map.dof_indices (test_elem.reference_elem, dof_indices);
97  uint n_dofs = uint(dof_indices.size());
98 
99  RealVectorX elem_solution = 5.3*RealVectorX::Ones(n_dofs);
100  test_elem.elem->set_solution(elem_solution);
101 
102  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(elem_solution),
103  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.elem->local_solution())));
104  }
105 
106  SECTION("Element shape can be transformed")
107  {
108  const Real V0 = test_elem.reference_elem->volume();
109 
110  // Stretch in x-direction
111  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0, 3.1, 1.0, 0.0, 0.0, 0.0);
112  REQUIRE(test_elem.reference_elem->volume() == 6.2);
113 
114  // Rotation about z-axis
115  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0, 1.0, 1.0, 0.0, 0.0, 60.0);
116  REQUIRE(test_elem.reference_elem->volume() == V0);
117 
118  // Rotation about y-axis
119  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0, 1.0, 1.0, 0.0, 30.0, 0.0);
120  REQUIRE(test_elem.reference_elem->volume() == V0);
121 
122  // Rotation about x-axis
123  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0, 1.0, 1.0, 20.0, 0.0, 0.0);
124  REQUIRE(test_elem.reference_elem->volume() == V0);
125 
126  // Shifted in x-direction
127  TEST::transform_element(test_elem.mesh, coords, 10.0, 0.0, 0.0, 1.0, 1.0, 0.0, 0.0, 0.0);
128  REQUIRE(test_elem.reference_elem->volume() == V0);
129 
130  // Shifted in y-direction
131  TEST::transform_element(test_elem.mesh, coords, 0.0, 7.5, 0.0, 1.0, 1.0, 0.0, 0.0, 0.0);
132  REQUIRE(test_elem.reference_elem->volume() == V0);
133 
134  // Shifted in z-direction
135  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 4.2, 1.0, 1.0, 0.0, 0.0, 0.0);
136  REQUIRE(test_elem.reference_elem->volume() == V0);
137  }
138 }
const MAST::NonlinearSystem & system() const
libMesh::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26
virtual unsigned int n_direct_strain_components()
row dimension of the direct strain matrix, also used for the bending operator row dimension ...
MAST::NonlinearImplicitAssembly assembly
libMesh::Elem * reference_elem
Pointer to the actual libMesh element object.
Definition: mast_mesh.h:51
libMesh::Real Real
virtual bool if_isotropic() const =0
return true if the property is isotropic
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...
TEST_CASE("structural_element_1d_base_tests", "[1D],[structural],[base]")
std::vector< double > eigen_matrix_to_std_vector(RealMatrixX M)
Converts an Eigen Matrix object to a std::vector.
virtual bool if_incompatible_modes() const
const MAST::ElementPropertyCardBase & elem_property()
returns a constant reference to the finite element object
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...
Matrix< Real, Dynamic, 1 > RealVectorX
const RealVectorX & local_solution(bool if_sens=false) const
virtual unsigned int n_von_karman_strain_components()
row dimension of the von Karman strain matrix
libMesh::ReplicatedMesh mesh
The actual libMesh mesh object.
Definition: mast_mesh.h:52