21 #include "libmesh/libmesh.h" 22 #include "libmesh/elem.h" 23 #include "libmesh/dof_map.h" 46 "[1D],[structural],[edge],[edge2],[linear]")
49 coords << -1.0, 1.0, 0.0,
57 double val_margin = (test_elem.
jacobian0.array().abs()).mean() * 1.490116119384766e-08;
59 libMesh::out <<
"J =\n" << test_elem.
jacobian0 << std::endl;
62 SECTION(
"Internal Jacobian (stiffness matrix) finite difference check")
67 val_margin = (test_elem.
jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
74 SECTION(
"Internal Jacobian (stiffness matrix) should be symmetric")
81 SECTION(
"Determinant of undeformed internal Jacobian (stiffness matrix) should be zero")
83 REQUIRE(test_elem.
jacobian0.determinant() == Approx(0.0).margin(1e-06));
87 SECTION(
"Internal Jacobian (stiffness matrix) eigenvalue check")
97 SelfAdjointEigenSolver<RealMatrixX> eigensolver(test_elem.
jacobian0,
false);
98 RealVectorX eigenvalues = eigensolver.eigenvalues();
99 libMesh::out <<
"Eigenvalues are:\n" << eigenvalues << std::endl;
101 for (uint i=0; i<eigenvalues.size(); i++)
103 if (std::abs(eigenvalues(i))<0.0001220703125)
113 REQUIRE(eigenvalues.minCoeff()>(-0.0001220703125));
117 SECTION(
"Internal Jacobian (stiffness matrix) is invariant to section orientation")
121 orientation(2) = 1.0;
131 SECTION(
"Internal Jacobian (stiffness matrix) is invariant to displacement solution")
134 elem_sol << 0.05727841, 0.08896581, 0.09541619, -0.03774913,
135 0.07510557, -0.07122266, -0.00979117, -0.08300009,
136 -0.03453369, -0.05487761, -0.01407677, -0.09268421;
145 SECTION(
"Internal Jacobian (stiffness matrix) is invariant to element x-location")
148 1.0, 1.0, 0.0, 0.0, 0.0);
157 SECTION(
"Internal Jacobian (stiffness matrix) is invariant to element y-location")
160 1.0, 1.0, 0.0, 0.0, 0.0);
169 SECTION(
"Internal Jacobian (stiffness matrix) is invariant to element z-location")
172 1.0, 1.0, 0.0, 0.0, 0.0);
181 SECTION(
"Internal Jacobian (stiffness matrix) checks for element aligned along y-axis")
189 new_coordinates << 0.0, 0.0, -1.0, 1.0, 0.0, 0.0;
196 val_margin = (test_elem.
jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
205 REQUIRE(test_elem.
jacobian.determinant() == Approx(0.0).margin(1e-06));
209 SECTION(
"Internal Jacobian (stiffness matrix) checks for element aligned along z-axis")
217 new_coordinates << 0.0, 0.0, 0.0, 0.0, -1.0, 1.0;
225 val_margin = (test_elem.
jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
234 REQUIRE(test_elem.
jacobian.determinant() == Approx(0.0).margin(1e-06));
238 SECTION(
"Internal Jacobian (stiffness matrix) checks for element rotated about z-axis")
242 1.0, 1.0, 0.0, 0.0, 63.4);
248 val_margin = (test_elem.
jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
257 REQUIRE(test_elem.
jacobian.determinant() == Approx(0.0).margin(1e-06));
261 SECTION(
"Internal Jacobian (stiffness matrix) checks for element rotated about y-axis")
265 1.0, 1.0, 0.0, 35.8, 0.0);
271 val_margin = (test_elem.
jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
280 REQUIRE(test_elem.
jacobian.determinant() == Approx(0.0).margin(1e-06));
284 SECTION(
"Internal Jacobian (stiffness matrix) checks for element stretched in x-direction")
287 3.2, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0);
293 val_margin = (test_elem.
jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
302 REQUIRE(test_elem.
jacobian.determinant() == Approx(0.0).margin(1e-06));
306 SECTION(
"Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, and rotated")
310 2.7, 6.4, 20.0, 47.8, -70.1);
316 val_margin = (test_elem.
jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
325 REQUIRE(test_elem.
jacobian.determinant() == Approx(0.0).margin(1e-06));
329 SECTION(
"Internal Jacobian (stiffness matrix) checks for element arbitrarily scaled, stretched, rotated, and displaced")
333 4.2, 1.5, -18.0, -24.8, 30.1);
337 elem_sol << 0.08158724, 0.07991906, -0.00719128, 0.02025461,
338 -0.04602193, 0.05280159, 0.03700081, 0.04636344,
339 0.05559377, 0.06448206, 0.08919238, -0.03079122;
347 val_margin = (test_elem.
jacobian_fd.array().abs()).mean() * 1.490116119384766e-08;
356 REQUIRE(test_elem.
jacobian.determinant() == Approx(0.0).margin(1e-06));
RealMatrixX jacobian
Matrix storage for Jacobian of the element in a perturbed/modified state.
virtual void init()
Only used by 1D sections.
RealVectorX elem_solution
libMesh::Elem * reference_elem
Pointer to the actual libMesh element object.
RealVectorX & y_vector()
returns value of the property val.
MAST::StructuralElement1D * elem
MAST::Solid1DSectionElementPropertyCard section
TEST_CASE("edge2_linear_structural", "[1D],[structural],[edge],[edge2],[linear]")
int n_nodes
Number of nodes per element in the test mesh.
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.
void approximate_internal_jacobian_with_finite_difference(MAST::StructuralElementBase &elem, const RealVectorX &initial_elem_solution, RealMatrixX &jacobian)
Approximates the internal Jacobian of an element using a 6th order accurate central finite difference...
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.
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 update_residual_and_jacobian0()
libMesh::LibMeshInit * p_global_init
void update_residual_and_jacobian()
void update_coordinates(RealMatrixX &new_coordinates)
Update the nodal coordinates in the mesh.
libMesh::ReplicatedMesh mesh
The actual libMesh mesh object.