35 #ifndef THREED_BEAM_FEA_H
36 #define THREED_BEAM_FEA_H
38 #ifdef EIGEN_USE_MKL_ALL
39 #include <Eigen/PardisoSupport>
42 #include <Eigen/SparseLU>
47 #include <Eigen/Geometry>
48 #include <Eigen/SparseCore>
70 typedef Eigen::Matrix<double, 12, 12, Eigen::RowMajor>
LocalMatrix;
110 SparseKelem.resize(12, 12);
111 SparseKelem.reserve(40);
125 void operator()(SparseMat &Kg,
const Job &job,
const std::vector<Tie> &ties);
143 void calcAelem(
const Eigen::Vector3d &nx,
const Eigen::Vector3d &nz);
170 SparseMat SparseKelem;
187 void loadBCs(SparseMat &Kg, ForceVector &force_vec,
const std::vector<BC> &BCs,
unsigned int num_nodes);
199 void loadTies(std::vector<Eigen::Triplet<double> > &triplets,
const std::vector<Tie> &ties);
210 std::vector<std::vector<double> >
computeTieForces(
const std::vector<Tie> &ties,
211 const std::vector<std::vector<double> > &nodal_displacements);
219 void loadForces(SparseMat &force_vec,
const std::vector<Force> &forces);
234 Summary
solve(
const Job &job,
235 const std::vector<BC> &BCs,
236 const std::vector<Force> &forces,
237 const std::vector<Tie> &ties,
238 const Options &options);
241 #endif // THREED_BEAM_FEA_H
Assembles the global stiffness matrix.
Definition: threed_beam_fea.h:97
GlobalStiffAssembler()
Default constructor.
Definition: threed_beam_fea.h:105
Eigen::Matrix< double, 12, 12, Eigen::RowMajor > LocalMatrix
Definition: threed_beam_fea.h:70
void calcAelem(const Eigen::Vector3d &nx, const Eigen::Vector3d &nz)
Updates the rotation and transposed rotation matrices.
Definition: threed_beam_fea.cpp:136
double norm(const Node &n1, const Node &n2)
Calculates the distance between 2 nodes.
Definition: threed_beam_fea.cpp:49
Summary solve(const Job &job, const std::vector< BC > &BCs, const std::vector< Force > &forces, const std::vector< Tie > &ties, const Options &options)
Solves the finite element analysis defined by the input Job, boundary conditions, and prescribed noda...
Definition: threed_beam_fea.cpp:378
LocalMatrix getKelem()
Returns the currently stored elemental stiffness matrix.
Definition: threed_beam_fea.h:149
void loadBCs(SparseMat &Kg, ForceVector &force_vec, const std::vector< BC > &BCs, unsigned int num_nodes)
Loads the boundary conditions into the global stiffness matrix and force vector.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > GlobalStiffMatrix
Definition: threed_beam_fea.h:65
void loadTies(std::vector< Eigen::Triplet< double > > &triplets, const std::vector< Tie > &ties)
Loads any tie constraints into the set of triplets that will become the global stiffness matrix...
Definition: threed_beam_fea.cpp:309
Eigen::SparseMatrix< double > SparseMat
Definition: threed_beam_fea.h:81
LocalMatrix getAelem()
Returns the currently stored rotation matrix.
Definition: threed_beam_fea.h:157
void calcKelem(unsigned int i, const Job &job)
Updates the elemental stiffness matrix for the ith element.
Definition: threed_beam_fea.cpp:54
void loadForces(SparseMat &force_vec, const std::vector< Force > &forces)
Loads the prescribed forces into the force vector.
Definition: threed_beam_fea.cpp:368
std::vector< std::vector< double > > computeTieForces(const std::vector< Tie > &ties, const std::vector< std::vector< double > > &nodal_displacements)
Computes the forces in the tie elements based on the nodal displacements of the FE analysis and the s...
Definition: threed_beam_fea.cpp:343
Eigen::Vector3d Node
A node that describes a mesh. Uses Eigen's predefined Vector class for added functionality.
Definition: containers.h:56
Contains a node list, element list, and the properties of each element.
Definition: containers.h:246
Eigen::Matrix< double, Eigen::Dynamic, 1 > ForceVector
Definition: threed_beam_fea.h:76
Definition: containers.h:41
void operator()(SparseMat &Kg, const Job &job, const std::vector< Tie > &ties)
Assembles the global stiffness matrix.
Definition: threed_beam_fea.cpp:226