#include <beamel2d.h>
Public Member Functions | |
void | beam_transf_matrix (long eid, vector &x, vector &z, matrix &tmat) |
beamel2d (void) | |
void | bf_matrix (matrix &n, double xi, double l, double kappa) |
void | dbf_matrix (matrix &n, double xi, double l, double kappa) |
void | define_meaning (long eid) |
void | geom_matrix (matrix &gm, double xi, double l, double kappa) |
void | initstr_matrix (long eid, long ri, long ci, matrix &ism) |
void | initstr_matrix_expl (long eid, long ri, long ci, matrix &ism) |
void | internal_forces (long lcid, long eid, vector &ifor) |
void | mass_matrix (long eid, long ri, long ci, matrix &mm) |
void | mass_matrix_expl (long eid, long ri, long ci, matrix &mm) |
void | nodal_displ (long lcid, long eid) |
void | nodal_forces (long lcid, long eid) |
void | nodeforces (long eid, long *le, double *nv, vector &nf) |
void | res_internal_forces (long lcid, long eid, vector &ifor) |
void | res_mass_matrix (long eid, matrix &mm) |
void | res_stiffness_matrix (long eid, matrix &sm) |
void | stiffness_matrix (long eid, long ri, long ci, matrix &sm) |
void | stiffness_matrix_expl (long eid, long ri, long ci, matrix &sm) |
void | stiffness_matrix_expl_local (long eid, long ri, long ci, matrix &sm) |
void | transf_matrix (ivector &nodes, matrix &tmat) |
~beamel2d (void) | |
Public Attributes | |
long * | cncomp |
array containing cumulative numbers of components of stress and strain tensors | |
long | intordism |
order of integration of initial stress matrix | |
long | intordmm |
order of integration of mass matrix | |
long ** | intordsm |
order of integration of stiffness matrix | |
long | napfun |
number of approximated functions on the element | |
long | nb |
number of blocks | |
long * | ncomp |
array containing numbers of components of stress and strain tensors | |
long | ndofe |
number of DOFs on the element | |
long ** | nip |
array of numbers of integration points in sets | |
long | nne |
number of nodes on one element | |
strastrestate | ssst |
stress/strain state | |
long | tncomp |
total number of components of stress and strain tensors | |
long | tnip |
class beamel2d defines beam element for 2D analysis
JK
Definition at line 16 of file beamel2d.h.
beamel2d | ( | void | ) |
~beamel2d | ( | void | ) |
function assembles transformation matrix from local element coordinate system to the global problem coordinate system x_g = T x_l
this beam element is formulated in x-z plane, x axis is identical with beam axis and it is oriented horizontally, z axis is oriented vertically (downwards)
x,z | - node coordinates | |
tmat | - transformation matrix |
JK, 3.2.2002
Definition at line 225 of file beamel2d.cpp.
References fillm(), and print_err().
Referenced by initstr_matrix(), initstr_matrix_expl(), mass_matrix(), mass_matrix_expl(), nodal_displ(), nodal_forces(), nodeforces(), stiffness_matrix(), and stiffness_matrix_expl().
void bf_matrix | ( | matrix & | n, | |
double | xi, | |||
double | l, | |||
double | kappa | |||
) |
function assembles matrix of approximation functions
ordering of approximation functions u - displacement along x w - deflection - displacement along z phi - rotation around y
n | - array containing matrix of approximation functions | |
xi | - natural coordinate from segment <0;1> | |
l | - length of the beam | |
kappa | - 6.0*E*I/k/G/A/l/l |
JK, 20.2.2002
Definition at line 88 of file beamel2d.cpp.
References vector::a, defl2D_fun(), dilat(), fillm(), and roty_fun().
Referenced by mass_matrix().
void dbf_matrix | ( | matrix & | n, | |
double | xi, | |||
double | l, | |||
double | kappa | |||
) |
function assembles matrix of derivatives of approximation function of deflection
n | - matrix of derivatives | |
xi | - natural coordinate from segment <0;1> | |
l | - length of the beam | |
kappa | - shear coefficient 6.0*E*I/k/G/A/l/l |
JK
Definition at line 122 of file beamel2d.cpp.
References vector::a, der_defl2D_fun(), and fillm().
Referenced by initstr_matrix().
void define_meaning | ( | long | eid | ) |
function defines meaning of DOFs at nodes
eid | - element id |
1.2.2005, JK
Definition at line 1076 of file beamel2d.cpp.
References mechtop::give_code_numbers(), mechtop::give_elemnodes(), node::meaning, Mt, ndofe, nne, nod, and mechtop::nodes.
Referenced by mechtop::define_meaning().
void geom_matrix | ( | matrix & | gm, | |
double | xi, | |||
double | l, | |||
double | kappa | |||
) |
function assembles strain-displacement (geometric) matrix
ordering of strain components du/dx - strain e_xx phi+dw/dx - strain e_xz dphi/dx - curvature
gm | - geometric matrix | |
xi | - natural coordinate from segment <0;1> | |
l | - length of the beam | |
kappa | - 6.0*E*I/k/G/A/l/l |
JK, 20.2.2002
Definition at line 152 of file beamel2d.cpp.
References vector::a, der_defl2D_fun(), der_dilat(), der_roty_fun(), fillm(), and roty_fun().
Referenced by stiffness_matrix().
void initstr_matrix | ( | long | eid, | |
long | ri, | |||
long | ci, | |||
matrix & | ism | |||
) |
function computes consistent geometric-stiffness matrix notation based on R. Clough, J. Penzien: Dynamic of Structures, McGraw-Hill, 2nd edition, 1993, pages 195-196
notation initial stress matrix is also used
eid | - number of element | |
ri,ci | - row and column indices | |
ism | - inital stress matrix |
JK
Definition at line 697 of file beamel2d.cpp.
References matrix::a, vector::a, beam_transf_matrix(), dbf_matrix(), mechtop::elements, fillm(), g, gauss_points(), mechcrsec::give_area(), mechtop::give_elemnodes(), mechcrsec::give_mominer(), mechtop::give_node_coord2dxz(), mechcrsec::give_shearcoeff(), glmatrixtransf(), intordism, element::ipp, lgmatrixtransf(), mechtop::locsystems(), mechmat::matstiff(), Mc, Mm, Mt, ndofe, nne, nnj(), nodes, print_err(), tncomp, and transf_matrix().
void initstr_matrix_expl | ( | long | eid, | |
long | ri, | |||
long | ci, | |||
matrix & | ism | |||
) |
function computes consistent geometric-stiffness matrix notation based on R. Clough, J. Penzien: Dynamic of Structures, McGraw-Hill, 2nd edition, 1993, pages 195-196
notation initial stress matrix is also used
initial stress matrix is formulated explicitly
eid | - number of element | |
ri,ci | - row and column indices | |
ism | - inital stress matrix |
JK
Definition at line 774 of file beamel2d.cpp.
References vector::a, beam_transf_matrix(), mechtop::elements, fillm(), g, gauss_points(), mechcrsec::give_area(), mechtop::give_elemnodes(), mechcrsec::give_mominer(), mechtop::give_node_coord2dxz(), mechcrsec::give_shearcoeff(), glmatrixtransf(), intordism, element::ipp, lgmatrixtransf(), mechtop::locsystems(), mechmat::matstiff(), Mc, Mm, Mt, ndofe, nne, nodes, print_err(), tncomp, and transf_matrix().
Referenced by initstiffmat().
void internal_forces | ( | long | lcid, | |
long | eid, | |||
vector & | ifor | |||
) |
Function computes internal forces (nodal forces and moments) that are expressed in global problem or local nodal coordinate system.
lcid | - load case id | |
eid | - element id | |
ri,ci | - row and column indices | |
ifor | - vector of internal forces |
JK, 12.8.2001
Definition at line 1013 of file beamel2d.cpp.
References addv(), eldispl(), fillv(), mxv(), ndofe, nne, nodes, and stiffness_matrix_expl().
Referenced by res_internal_forces().
void mass_matrix | ( | long | eid, | |
long | ri, | |||
long | ci, | |||
matrix & | mm | |||
) |
function computes consistent mass matrix of 2D beam element influence of inertial forces from rotations is accounted
eid | - number of element | |
ri,ci | - row and column indices | |
mm | - mass matrix |
JK, 3.2.2006
Definition at line 498 of file beamel2d.cpp.
References matrix::a, vector::a, bdbj(), beam_transf_matrix(), bf_matrix(), mechtop::elements, fillm(), g, gauss_points(), mechcrsec::give_area(), mechcrsec::give_densitye(), mechtop::give_elemnodes(), mechcrsec::give_mominer(), mechtop::give_node_coord2dxz(), mechcrsec::give_shearcoeff(), glmatrixtransf(), intordmm, element::ipp, lgmatrixtransf(), mechtop::locsystems(), mechmat::matstiff(), Mc, Mm, Mt, napfun, ndofe, nne, nodes, print_err(), tncomp, and transf_matrix().
void mass_matrix_expl | ( | long | eid, | |
long | ri, | |||
long | ci, | |||
matrix & | mm | |||
) |
function computes consistent mass matrix of 2D beam element influence of inertial forces from rotations is accounted
eid | - number of element | |
ri,ci | - row and column indices | |
mm | - mass matrix |
JK, 3.2.2006
Definition at line 590 of file beamel2d.cpp.
References beam_transf_matrix(), mechtop::elements, fillm(), g, mechcrsec::give_area(), mechcrsec::give_densitye(), mechtop::give_elemnodes(), mechcrsec::give_mominer(), mechtop::give_node_coord2dxz(), mechcrsec::give_shearcoeff(), glmatrixtransf(), element::ipp, lgmatrixtransf(), mechtop::locsystems(), mechmat::matstiff(), Mc, Mm, Mt, ndofe, nne, nodes, print_err(), tncomp, and transf_matrix().
Referenced by nodal_forces(), and res_mass_matrix().
void nodal_displ | ( | long | lcid, | |
long | eid | |||
) |
function stores end displacements and rotations to integration points
lcid | - load case id | |
eid | - element id |
JK, 24.5.2006
Definition at line 861 of file beamel2d.cpp.
References beam_transf_matrix(), copyv(), eldispl(), mechtop::elements, f, mechtop::give_elemnodes(), mechtop::give_node_coord2dxz(), glvectortransf(), element::ipp, lgvectortransf(), mechtop::locsystems(), Mm, Mt, ndofe, nne, nodes, mechmat::storestrain(), and transf_matrix().
Referenced by compute_ipstrains(), and compute_nodestrains().
void nodal_forces | ( | long | lcid, | |
long | eid | |||
) |
function computes nodal forces and moments expressed in local coordinate system this is equivalent to functions computing stresses nodal forces and moments are better quantities in case of beams
lcid | - load case id | |
eid | - element id |
JK, 20.2.2002, revised 2.9.2006
Definition at line 918 of file beamel2d.cpp.
References beam_transf_matrix(), cmulv(), copyv(), loadel::eid, eigen_dynamics, eldispl(), mechtop::elements, f, forced_dynamics, mechtop::give_elemnodes(), mechtop::give_node_coord2dxz(), glvectortransf(), element::ipp, mechbclc::lc, lgvectortransf(), mechtop::locsystems(), loadcase::loe, Lsrs, mass_matrix_expl(), Mb, Mm, Mp, Mt, mxv(), ndofe, loadel::nf, loadcase::nle, nne, nodes, stiffness_matrix_expl(), mechmat::storestress(), probdesc::tprob, transf_matrix(), and lhsrhs::w.
Referenced by compute_ipstresses(), compute_nodestresses(), and computestresses().
void nodeforces | ( | long | eid, | |
long * | le, | |||
double * | nv, | |||
vector & | nf | |||
) |
function assembles end forces and moments cuased by a load acting on element
eid | - element id | |
le | - type of coordinate system | |
nv | - values defining load | |
nf | - vector of end forces and moments |
30. 9. 2012
Definition at line 1104 of file beamel2d.cpp.
References ivector::a, beam_transf_matrix(), gelement::cne, condense_vector(), mechtop::elements, gtopology::gelements, mechcrsec::give_area(), gtopology::give_cn(), mechcrsec::give_mominer(), mechtop::give_node_coord2dxz(), Gtm, h(), element::ipp, lgvectortransfblock(), mechmat::matstiff(), Mc, Mm, Mt, ndofe, nne, stiffness_matrix_expl(), and tncomp.
Referenced by loadel::edgeload().
void res_internal_forces | ( | long | lcid, | |
long | eid, | |||
vector & | ifor | |||
) |
function computes internal forces
lcid | - load case id | |
eid | - element id | |
ifor | - vector of internal forces |
JK, 12.8.2001
Definition at line 1064 of file beamel2d.cpp.
References internal_forces().
Referenced by elem_internal_forces(), and elem_nonloc_internal_forces().
void res_mass_matrix | ( | long | eid, | |
matrix & | mm | |||
) |
function computes consistent mass matrix of 2D beam element influence of inertial forces from rotations is accounted
eid | - number of element | |
mm | - mass matrix |
JK, 3.2.2006
Definition at line 574 of file beamel2d.cpp.
References mass_matrix_expl().
Referenced by massmat().
void res_stiffness_matrix | ( | long | eid, | |
matrix & | sm | |||
) |
function computes stiffness matrix of two-dimensional beam finite element element is based on Mindlin theory (shear stresses are taken into account)
eid | - number of element | |
sm | - stiffness matrix |
JK, 3.2.2002
Definition at line 337 of file beamel2d.cpp.
References stiffness_matrix_expl().
Referenced by stiffmat().
void stiffness_matrix | ( | long | eid, | |
long | ri, | |||
long | ci, | |||
matrix & | sm | |||
) |
function computes stiffness matrix of two-dimensional beam finite element element is based on Mindlin theory (shear stresses are taken into account)
eid | - number of element | |
ri,ci | - row and column indices | |
sm | - stiffness matrix |
JK, 3.2.2002
Definition at line 263 of file beamel2d.cpp.
References matrix::a, vector::a, bdbj(), beam_transf_matrix(), mechtop::elements, fillm(), g, gauss_points(), geom_matrix(), mechcrsec::give_area(), mechtop::give_elemnodes(), mechcrsec::give_mominer(), mechtop::give_node_coord2dxz(), mechcrsec::give_shearcoeff(), glmatrixtransf(), intordsm, element::ipp, lgmatrixtransf(), mechtop::locsystems(), mechmat::matstiff(), Mc, Mm, Mt, ndofe, nne, nodes, print_err(), tncomp, and transf_matrix().
void stiffness_matrix_expl | ( | long | eid, | |
long | ri, | |||
long | ci, | |||
matrix & | sm | |||
) |
Definition at line 422 of file beamel2d.cpp.
References ivector::a, beam_transf_matrix(), gelement::cne, condense_matrix(), mechtop::elements, fillm(), g, gtopology::gelements, mechcrsec::give_area(), gtopology::give_cn(), mechtop::give_elemnodes(), mechcrsec::give_mominer(), mechtop::give_node_coord2dxz(), mechcrsec::give_shearcoeff(), glmatrixtransf(), Gtm, element::ipp, lgmatrixtransf(), mechtop::locsystems(), mechmat::matstiff(), Mc, Mm, Mt, ndofe, nne, nodes, stiffness_matrix_expl_local(), tncomp, and transf_matrix().
Referenced by internal_forces(), nodal_forces(), nodeforces(), and res_stiffness_matrix().
void stiffness_matrix_expl_local | ( | long | eid, | |
long | ri, | |||
long | ci, | |||
matrix & | sm | |||
) |
function computes stiffness matrix of two-dimensional beam finite element element is based on Mindlin theory (shear stresses are taken into account) stiffness matrix is formulated explicitly
eid | - number of element | |
ri,ci | - row and column indices | |
sm | - stiffness matrix |
JK, 3.2.2002
Definition at line 354 of file beamel2d.cpp.
References mechtop::elements, fillm(), g, mechcrsec::give_area(), mechcrsec::give_mominer(), mechtop::give_node_coord2dxz(), mechcrsec::give_shearcoeff(), element::ipp, mechmat::matstiff(), Mc, Mm, Mt, ndofe, nne, nodes, print_err(), and tncomp.
Referenced by stiffness_matrix_expl().
function assembles transformation matrix from local nodal coordinate system to the global coordinate system x_g = T x_l
nodes | - array containing node numbers | |
tmat | - transformation matrix |
JK, 3.2.2002
Definition at line 190 of file beamel2d.cpp.
References node::e1, node::e2, fillm(), matrix::m, Mt, ivector::n, mechtop::nodes, and node::transf.
Referenced by initstr_matrix(), initstr_matrix_expl(), mass_matrix(), mass_matrix_expl(), nodal_displ(), nodal_forces(), stiffness_matrix(), and stiffness_matrix_expl().
long* cncomp |
array containing cumulative numbers of components of stress and strain tensors
Definition at line 57 of file beamel2d.h.
Referenced by beamel2d(), and ~beamel2d().
long intordism |
order of integration of initial stress matrix
Definition at line 65 of file beamel2d.h.
Referenced by beamel2d(), initstr_matrix(), and initstr_matrix_expl().
long intordmm |
order of integration of mass matrix
Definition at line 63 of file beamel2d.h.
Referenced by beamel2d(), and mass_matrix().
long** intordsm |
order of integration of stiffness matrix
Definition at line 61 of file beamel2d.h.
Referenced by beamel2d(), mechtop::give_intordsm(), stiffness_matrix(), and ~beamel2d().
long napfun |
number of approximated functions on the element
Definition at line 59 of file beamel2d.h.
Referenced by beamel2d(), mechtop::give_napfun(), and mass_matrix().
long nb |
number of blocks
Definition at line 69 of file beamel2d.h.
Referenced by beamel2d(), mechtop::give_nb(), mechtop::give_nb_te(), and ~beamel2d().
long* ncomp |
array containing numbers of components of stress and strain tensors
Definition at line 55 of file beamel2d.h.
Referenced by beamel2d(), and ~beamel2d().
long ndofe |
number of DOFs on the element
Definition at line 49 of file beamel2d.h.
Referenced by beamel2d(), define_meaning(), mechtop::give_ndofe(), initstr_matrix(), initstr_matrix_expl(), internal_forces(), mass_matrix(), mass_matrix_expl(), nodal_displ(), nodal_forces(), nodeforces(), stiffness_matrix(), stiffness_matrix_expl(), and stiffness_matrix_expl_local().
long** nip |
array of numbers of integration points in sets
Definition at line 67 of file beamel2d.h.
Referenced by beamel2d(), mechtop::give_nip(), and ~beamel2d().
long nne |
number of nodes on one element
Definition at line 51 of file beamel2d.h.
Referenced by beamel2d(), define_meaning(), mechtop::give_nne(), initstr_matrix(), initstr_matrix_expl(), internal_forces(), mass_matrix(), mass_matrix_expl(), nodal_displ(), nodal_forces(), nodeforces(), stiffness_matrix(), stiffness_matrix_expl(), and stiffness_matrix_expl_local().
stress/strain state
Definition at line 71 of file beamel2d.h.
Referenced by beamel2d(), and mechtop::give_ssst().
long tncomp |
total number of components of stress and strain tensors
Definition at line 53 of file beamel2d.h.
Referenced by beamel2d(), mechtop::give_ncomp(), mechtop::give_tncomp(), initstr_matrix(), initstr_matrix_expl(), mass_matrix(), mass_matrix_expl(), nodeforces(), stiffness_matrix(), stiffness_matrix_expl(), and stiffness_matrix_expl_local().
long tnip |
Definition at line 73 of file beamel2d.h.
Referenced by beamel2d(), and mechtop::give_tnip().