MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mast_quad4_linear_structural_strain_displacement_matrix.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 #define protected public
21 
22 // libMesh includes
23 #include "libmesh/libmesh.h"
24 #include "libmesh/elem.h"
25 #include "libmesh/equation_systems.h"
26 #include "libmesh/dof_map.h"
27 
28 // MAST includes
29 #include "base/parameter.h"
37 #include "base/nonlinear_system.h"
38 #include "mesh/geom_elem.h"
39 #include "mesh/fe_base.h"
41 
42 // Test includes
43 #include "catch.hpp"
44 #include "test_helpers.h"
46 
47 extern libMesh::LibMeshInit* p_global_init;
48 
55 TEST_CASE("quad4_linear_structural_strain_displacement_matrix",
56  "[quad],[quad4],[linear],][structural],[2D],[element]")
57 {
58  RealMatrixX coords = RealMatrixX::Zero(3,4);
59  coords << -1.0, 1.0, 1.0, -1.0,
60  -1.0, -1.0, 1.0, 1.0,
61  0.0, 0.0, 0.0, 0.0;
62  TEST::TestStructuralSingleElement2D test_elem(libMesh::QUAD4, coords);
63 
64  const Real V0 = test_elem.reference_elem->volume();
65  REQUIRE(test_elem.reference_elem->volume() == TEST::get_shoelace_area(coords));
66 
67  // Set the strain type to linear for the section
69 
70 
71  std::unique_ptr<MAST::FEBase> fe(test_elem.geom_elem.init_fe(true, false,
72  test_elem.section.extra_quadrature_order(test_elem.geom_elem)));
73 
74  const uint n_dofs_per_node = 6;
76  Bmat_lin,
77  Bmat_nl_x,
78  Bmat_nl_y,
79  Bmat_nl_u,
80  Bmat_nl_v,
81  Bmat_bend,
82  Bmat_vk;
83 
84  const uint n_phi = (unsigned int)fe->get_phi().size();
85  const uint n1 = test_elem.elem->n_direct_strain_components();
86  const uint n2 = 6*n_phi;
87  const uint n3 = test_elem.elem->n_von_karman_strain_components();
88 
89  RealMatrixX mat_x = RealMatrixX::Zero(3,2);
90  RealMatrixX mat_y = RealMatrixX::Zero(3,2);
91 
92  RealVectorX strain = RealVectorX::Zero(3);
93 
94  uint qp = 0;
95 
100  Bmat_lin.reinit(n1, test_elem.structural_system.n_vars(), n_phi); // three stress-strain components
101  Bmat_nl_x.reinit(2, test_elem.structural_system.n_vars(), n_phi);
102  Bmat_nl_y.reinit(2, test_elem.structural_system.n_vars(), n_phi);
103  Bmat_nl_u.reinit(2, test_elem.structural_system.n_vars(), n_phi);
104  Bmat_nl_v.reinit(2, test_elem.structural_system.n_vars(), n_phi);
105 
106  test_elem.elem->initialize_green_lagrange_strain_operator(qp, *fe, test_elem.elem_solution,
107  strain, mat_x, mat_y, Bmat_lin, Bmat_nl_x, Bmat_nl_y, Bmat_nl_u,
108  Bmat_nl_v);
109 
115  Bmat_bend.reinit(n1, test_elem.structural_system.n_vars(), n_phi);
116 
122  Bmat_vk.reinit(n3, test_elem.structural_system.n_vars(), n_phi); // only dw/dx and dw/dy
123 
124 
125  const std::vector<std::vector<libMesh::RealVectorValue>>& dphi = fe->get_dphi();
126 
127  SECTION("Linear in-plane strain displacement matrix size")
128  {
129  REQUIRE( Bmat_lin.m() == 3 ); // three strains, e_xx, e_yy, e_xy
130  REQUIRE( Bmat_lin.n() == test_elem.n_nodes * n_dofs_per_node);
131  }
132 
133  SECTION("Linear in-plane strain displacement matrix values")
134  {
135  // TODO: This requires the vector_mult method be working correctly. Is that ok?
136 
137  // First get a RealMatrixX representation of this matrix
138  uint m = Bmat_lin.m();
139  uint n = Bmat_lin.n();
140  RealMatrixX Bmat_lin_mat = RealMatrixX::Zero(m,n);
141  for (uint i=0; i<n; i++)
142  {
143  RealVectorX Ivec = RealVectorX::Zero(n);
144  RealVectorX result = RealVectorX::Zero(m);
145  Ivec(i) = 1.0;
146  Bmat_lin.vector_mult(result, Ivec);
147  Bmat_lin_mat.col(i) = result;
148  }
149 
158  RealMatrixX Bmat_lin_true = RealMatrixX::Zero(m,n);
159  for (uint i=0; i<test_elem.n_nodes; i++)
160  {
161  Bmat_lin_true(0,i) = dphi[i][qp](0);
162  Bmat_lin_true(2,i+test_elem.n_nodes) = dphi[i][qp](0);
163 
164  Bmat_lin_true(1,i+test_elem.n_nodes) = dphi[i][qp](1);
165  Bmat_lin_true(2,i) = dphi[i][qp](1);
166  }
167 
168  std::vector<double> Bmat_lin_test = TEST::eigen_matrix_to_std_vector(Bmat_lin_mat);
169  std::vector<double> Bmat_lin_required = TEST::eigen_matrix_to_std_vector(Bmat_lin_true);
170 
171  REQUIRE_THAT(Bmat_lin_test, Catch::Approx<double>(Bmat_lin_required));
172  }
173 }
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
void reinit(unsigned int n_interpolated_vars, unsigned int n_discrete_vars, unsigned int n_discrete_dofs_per_var)
this initializes the operator for number of rows and variables, assuming that all variables has the s...
MAST::StructuralSystemInitialization structural_system
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
int n_nodes
Number of nodes per element in the test mesh.
Definition: mast_mesh.h:49
unsigned int n() const
virtual std::unique_ptr< MAST::FEBase > init_fe(bool init_grads, bool init_second_order_derivative, int extra_quadrature_order=0) const
initializes the finite element shape function and quadrature object with the order of quadrature rule...
Definition: geom_elem.cpp:148
std::vector< double > eigen_matrix_to_std_vector(RealMatrixX M)
Converts an Eigen Matrix object to a std::vector.
libMesh::LibMeshInit * p_global_init
Definition: init_catch2.cpp:26
TEST_CASE("quad4_linear_structural_strain_displacement_matrix", "[quad],[quad4],[linear],][structural],[2D],[element]")
References https://studiumbook.com/properties-of-shape-function-fea/ https://www.ccg.msm.cam.ac.uk/images/FEMOR_Lecture_2.pdf
Matrix< Real, Dynamic, Dynamic > RealMatrixX
Matrix< Real, Dynamic, 1 > RealVectorX
void vector_mult(T &res, const T &v) const
res = [this] * v
unsigned int m() const
virtual void initialize_green_lagrange_strain_operator(const unsigned int qp, const MAST::FEBase &fe, const RealVectorX &local_disp, RealVectorX &epsilon, RealMatrixX &mat_x, RealMatrixX &mat_y, MAST::FEMOperatorMatrix &Bmat_lin, MAST::FEMOperatorMatrix &Bmat_nl_x, MAST::FEMOperatorMatrix &Bmat_nl_y, MAST::FEMOperatorMatrix &Bmat_nl_u, MAST::FEMOperatorMatrix &Bmat_nl_v)
virtual int extra_quadrature_order(const MAST::GeomElem &elem) const
returns the extra quadrature order (on top of the system) that this element should use...
virtual unsigned int n_von_karman_strain_components()
row dimension of the von Karman strain matrix