MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_structural_element_2d.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_2d_base_tests",
45  "[2D][structural][base]")
46 {
47  RealMatrixX coords = RealMatrixX::Zero(3,4);
48  coords << -1.0, 1.0, 1.0, -1.0,
49  -1.0, -1.0, 1.0, 1.0,
50  0.0, 0.0, 0.0, 0.0;
51  TEST::TestStructuralSingleElement2D test_elem(libMesh::QUAD4, coords);
52 
53  REQUIRE(test_elem.reference_elem->volume() == TEST::get_shoelace_area(coords));
54 
55  SECTION("Element returns proper number of strain components")
56  {
57  REQUIRE(test_elem.elem->n_direct_strain_components() == 3);
58  REQUIRE(test_elem.elem->n_von_karman_strain_components() == 2);
59  }
60 
61  SECTION("Incompatible modes flag returns false for basic element")
62  {
63  REQUIRE_FALSE(test_elem.elem->if_incompatible_modes());
64  }
65 
66  SECTION("Isotropic flag returns true for basic element")
67  {
68  const MAST::ElementPropertyCardBase& elem_section = test_elem.elem->elem_property();
69  CHECK( elem_section.if_isotropic() );
70  }
71 
72  SECTION("Check setting/getting element local solution")
73  {
74  const libMesh::DofMap& dof_map = test_elem.assembly.system().get_dof_map();
75  std::vector<libMesh::dof_id_type> dof_indices;
76  dof_map.dof_indices (test_elem.reference_elem, dof_indices);
77  uint n_dofs = uint(dof_indices.size());
78 
79  RealVectorX elem_solution = 5.3*RealVectorX::Ones(n_dofs);
80  test_elem.elem->set_solution(elem_solution);
81 
82  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(elem_solution),
83  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.elem->local_solution())));
84  }
85 
86  SECTION("Check setting/getting element local sensitivity solution")
87  {
88  const libMesh::DofMap& dof_map = test_elem.assembly.system().get_dof_map();;
89  std::vector<libMesh::dof_id_type> dof_indices;
90  dof_map.dof_indices (test_elem.reference_elem, dof_indices);
91  uint n_dofs = uint(dof_indices.size());
92 
93  RealVectorX elem_solution_sens = 3.1*RealVectorX::Ones(n_dofs);
94  test_elem.elem->set_solution(elem_solution_sens, true);
95 
96  const RealVectorX& local_solution_sens = test_elem.elem->local_solution(true);
97 
98  REQUIRE_THAT(TEST::eigen_matrix_to_std_vector(elem_solution_sens),
99  Catch::Approx<double>(TEST::eigen_matrix_to_std_vector(test_elem.elem->local_solution(true))));
100  }
101 
102  SECTION("Element shape can be transformed")
103  {
104  const Real V0 = test_elem.reference_elem->volume();
105 
106  // Stretch in x-direction
107  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0, 3.1, 1.0, 0.0, 0.0, 0.0);
108  REQUIRE(test_elem.reference_elem->volume() == 12.4);
109 
110  // Stretch in y-direction
111  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 0.0, 1.0, 3.1, 0.0, 0.0, 0.0);
112  REQUIRE(test_elem.reference_elem->volume() == 12.4);
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  // Shear in x
139  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 4.2, 1.0, 1.0, 0.0, 0.0, 0.0, 5.2, 0.0);
140  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
141 
142  // Shear in y
143  TEST::transform_element(test_elem.mesh, coords, 0.0, 0.0, 4.2, 1.0, 1.0, 0.0, 0.0, 0.0, 0.0, -6.4);
144  REQUIRE(test_elem.reference_elem->volume() == Approx(V0));
145  }
146 }
libMesh::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26
const MAST::NonlinearSystem & system() const
virtual unsigned int n_direct_strain_components()
row dimension of the direct strain matrix, also used for the bending operator row dimension ...
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.
libMesh::Real Real
virtual bool if_isotropic() const =0
return true if the property is isotropic
MAST::NonlinearImplicitAssembly assembly
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.
const MAST::ElementPropertyCardBase & elem_property()
returns a constant reference to the finite element object
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual bool if_incompatible_modes() const
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
TEST_CASE("structural_element_2d_base_tests", "[2D][structural][base]")
libMesh::ReplicatedMesh mesh
The actual libMesh mesh object.
Definition: mast_mesh.h:52
virtual unsigned int n_von_karman_strain_components()
row dimension of the von Karman strain matrix