beamel3d Class Reference

#include <beamel3d.h>

List of all members.

Public Member Functions

void beam_transf_matrix (matrix &tmat, double &dl, vector &vec, vector &x, vector &y, vector &z, long eid)
 beamel3d (void)
void bf_matrix (matrix &n, double s, double dl, double gy, double gz)
void bf_matrix_transl (matrix &n, double xi, double l, double gy, double gz)
void define_meaning (long eid)
void geom_matrix (matrix &n, double s, double dl, double gy, double gz)
void initstr_matrix (long eid, long ri, long ci, matrix &ism)
void internal_forces (long lcid, long eid, long ri, long ci, vector &ifor)
void load_matrix (long eid, matrix &lm)
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_local (long eid, long ri, long ci, matrix &sm)
void stiffness_matrix_local2 (long eid, long ri, long ci, matrix &sm)
void transf_matrix (ivector &nodes, matrix &tmat)
 ~beamel3d (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
 total number of integration points

Detailed Description

class beamel3d

PF

Definition at line 16 of file beamel3d.h.


Constructor & Destructor Documentation

beamel3d ( void   ) 

Definition at line 11 of file beamel3d.cpp.

References cncomp, intordism, intordmm, intordsm, napfun, nb, ncomp, ndofe, nip, nne, spacebeam, ssst, tncomp, and tnip.

~beamel3d ( void   ) 

Definition at line 63 of file beamel3d.cpp.

References cncomp, intordsm, nb, ncomp, and nip.


Member Function Documentation

void beam_transf_matrix ( matrix tmat,
double &  dl,
vector vec,
vector x,
vector y,
vector z,
long  eid 
)

The function assembles transformation matrix from local element coordinate system to global coordinate system x_g = T x_l.

Columns of the transformation matrix are created by coordinates of local base vectors expressed in global coordinate system.

Parameters:
tmat - transformation matrix
dl - lenght of the beam
vec - direction vector of local axis z (in global coordinates)
x,y,z - vectors of nodal coordinates in global system
eid - element id

JK, 3.2.2002, revised 1.9.2006

Definition at line 312 of file beamel3d.cpp.

References fillm(), and print_err().

Referenced by initstr_matrix(), internal_forces(), load_matrix(), mass_matrix(), nodal_displ(), nodal_forces(), stiffness_matrix(), stiffness_matrix_local(), and stiffness_matrix_local2().

void bf_matrix ( matrix n,
double  s,
double  dl,
double  gy,
double  gz 
)

function assembles matrix of approximation functions

ordering of approximation functions u v w - displacement along x,y,z fix,fiy,fiz - rotation around x,y,z

Parameters:
n - matrix of approximation functions
s - natural coordinate from segment <0;1>
dl - length of the beam
gy - 6.0*E*I_y/k/G/A/l/l
gz - 6.0*E*I_z/k/G/A/l/l

PF, 20.12.2002, revised by JK 2.9.2006

Definition at line 160 of file beamel3d.cpp.

References ll.

Referenced by load_matrix().

void bf_matrix_transl ( matrix n,
double  xi,
double  l,
double  gy,
double  gz 
)

function assembles matrix of approximation functions for translational DOFs

ordering of approximation functions u v w - displacement along x,y,z fix,fiy,fiz - rotation around x,y,z are neglected

ordering of nodal values u_1, v_1, w_1, phi_x1, phi_y1, phi_z1, u_2, v_2, w_2, phi_x2, phi_y2, phi_z2

Parameters:
n - matrix of approximation functions
xi - natural coordinate from segment <0;1>
l - length of the beam
gy - 6.0*E*I_y/k/G/A/l/l
gz - 6.0*E*I_z/k/G/A/l/l

JK, 16.11.2006

Definition at line 219 of file beamel3d.cpp.

References fillm().

Referenced by mass_matrix().

void define_meaning ( long  eid  ) 

The function defines meaning of DOFs at nodes.

Parameters:
eid - element id

1.2.2005, JK

Definition at line 1267 of file beamel3d.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 n,
double  s,
double  dl,
double  gy,
double  gz 
)

function assembles geometric matrix (function assembles matrix of derivatives of approximation functions)

ordering of approximation functions

du/dx - strain e_xx phi+dv/dx - strain e_yz phi+dw/dx - strain e_xz dphix/dx - curvature dphiy/dx - curvature dphiz/dx - curvature

Parameters:
n - array containing geometric matrix
s - natural coordinate from segment <0;1>
dl - length of the beam
gy - 6.0*E*I_y/k/G/A/l/l
gz - 6.0*E*I_z/k/G/A/l/l

PF, 20.12.2002, revised by JK 2.9.2006

Definition at line 101 of file beamel3d.cpp.

References fillm(), and ll.

Referenced by internal_forces().

void initstr_matrix ( long  eid,
long  ri,
long  ci,
matrix ism 
)
void internal_forces ( long  lcid,
long  eid,
long  ri,
long  ci,
vector ifor 
)
void load_matrix ( long  eid,
matrix lm 
)
void mass_matrix ( long  eid,
long  ri,
long  ci,
matrix mm 
)

function computes consistent mass matrix of 3D beam element rotational DOFs are neglected

Parameters:
eid - number of element
ri,ci - row and column indices
mm - mass matrix of one element

JK, 3.2.200

Definition at line 727 of file beamel3d.cpp.

References matrix::a, beam_transf_matrix(), bf_matrix_transl(), fillm(), gauss_points(), mechcrsec::give_area(), mechcrsec::give_densitye(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), mechcrsec::give_shearcoeff(), mechcrsec::give_vectorlcs(), glmatrixtransf(), intordmm, lgmatrixtransf(), mechtop::locsystems(), Mc, Mt, ndofe, nne, nnj(), nodes, and transf_matrix().

Referenced by res_mass_matrix().

void mass_matrix_expl ( long  eid,
long  ri,
long  ci,
matrix mm 
)
void nodal_displ ( long  lcid,
long  eid 
)

The function stores nodal displacements and rotations in the local element coordinate system. This is equivalent to functions computing strains. Nodal displacements and rotations are better quantities in case of beams.

Parameters:
eid - element id
lcid - load case id

JK, 20.2.2002

Definition at line 994 of file beamel3d.cpp.

References beam_transf_matrix(), copyv(), eldispl(), mechtop::elements, f, mechtop::give_elemnodes(), mechtop::give_node_coord3d(), mechcrsec::give_vectorlcs(), glvectortransf(), element::ipp, lgvectortransf(), mechtop::locsystems(), Mc, Mm, Mt, ndofe, nne, nodes, mechmat::storestrain(), and transf_matrix().

Referenced by compute_ipstrains(), and compute_nodestrains().

void nodal_forces ( long  lcid,
long  eid 
)

The function computes nodal forces and moments expressed in local coodinate system. This is equivalent to functions computing stresses on plane or space elements. Nodal forces and moments are better quantities in case of beams.

x-axis goes through the beam z-axis is defined in cross section y-axis is defined by right hand orientation

Parameters:
eid - element id
lcid - load case id

JK, 20.2.2002, revised 2.9.2006

Definition at line 1063 of file beamel3d.cpp.

References beam_transf_matrix(), copyv(), eldispl(), mechtop::elements, f, mechtop::give_elemnodes(), mechtop::give_node_coord3d(), mechcrsec::give_vectorlcs(), glvectortransf(), element::ipp, lgvectortransf(), mechtop::locsystems(), Mc, Mm, Mt, mxv(), ndofe, nne, nodes, stiffness_matrix(), mechmat::storestress(), and transf_matrix().

Referenced by compute_ipstresses(), compute_nodestresses(), and computestresses().

void nodeforces ( long  eid,
long *  le,
double *  nv,
vector nf 
)

The function calculates member forces due to continuous load of element.

Parameters:
eid - element id
le -
nv - nodal values of continuous load
nf - resulting vector of member forces
Returns:
The function returns resulting member forces in the parameter nf.

Created by Jan Zitny 10.2012

Definition at line 666 of file beamel3d.cpp.

References ivector::a, gelement::cne, condense_vector(), gtopology::gelements, gtopology::give_cn(), mechtop::give_node_coord3d(), Gtm, Mt, ndofe, nne, print_err(), and stiffness_matrix_local().

Referenced by loadel::edgeload().

void res_internal_forces ( long  lcid,
long  eid,
vector ifor 
)

The function computes internal forces, i.e. member forces of the element expressed in the global coordinate system or local nodal coordinate systems.

Parameters:
lcid - load case id
eid - element id
ri,ci - row and column indices
ifor - vector of internal forces

20.12.2002

Definition at line 1136 of file beamel3d.cpp.

References vector::a, eldispl(), mxv(), ndofe, and stiffness_matrix().

Referenced by elem_internal_forces(), and elem_nonloc_internal_forces().

void res_mass_matrix ( long  eid,
matrix mm 
)

function computes consistent mass matrix of 3D beam element

Parameters:
eid - element id
mm - mass matrix

JK, 2.9.2006

Definition at line 795 of file beamel3d.cpp.

References mass_matrix().

Referenced by massmat().

void res_stiffness_matrix ( long  eid,
matrix sm 
)

function computes stiffness matrix of 3D beam element

Parameters:
eid - number of element
sm - stiffness matrix

JK, 3.2.2002

Definition at line 432 of file beamel3d.cpp.

References stiffness_matrix().

Referenced by stiffmat().

void stiffness_matrix ( long  eid,
long  ri,
long  ci,
matrix sm 
)

The function computes stiffness matrix of 3D beam element expressed in the global coordinate system or local nodal coordinate systems. Additionally, the condensation of the stiffness matrix is performed.

Parameters:
eid - number of element
ri,ci - row and column indices
sm - stiffness matrix
Returns:
The function returns resulting matrix in the parameter sm.

PF, 20.12.2002 Jan Zitny 2012

Definition at line 608 of file beamel3d.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_coord3d(), mechcrsec::give_shearcoeff(), mechcrsec::give_vectorlcs(), glmatrixtransf(), Gtm, element::ipp, lgmatrixtransf(), ll, mechtop::locsystems(), mechmat::matstiff(), Mc, Mm, Mt, ndofe, nne, nodes, stiffness_matrix_local2(), tncomp, and transf_matrix().

Referenced by nodal_forces(), res_internal_forces(), and res_stiffness_matrix().

void stiffness_matrix_local ( long  eid,
long  ri,
long  ci,
matrix sm 
)

The function computes stiffness matrix of 3D beam element expressed in the local element coordinate system.

Parameters:
eid - number of element
ri,ci - row and column indices
sm - stiffness matrix
Returns:
The function returns resulting matrix in the parameter sm.

PF, 20.12.2002 Jan Zitny, 2012

Definition at line 453 of file beamel3d.cpp.

References beam_transf_matrix(), mechtop::elements, fillm(), g, mechcrsec::give_area(), mechtop::give_elemnodes(), mechcrsec::give_mominer(), mechtop::give_node_coord3d(), mechcrsec::give_shearcoeff(), mechcrsec::give_vectorlcs(), element::ipp, ll, mechmat::matstiff(), Mc, Mm, Mp, Mt, ndofe, nne, nodes, tncomp, and probdesc::zero.

Referenced by nodeforces().

void stiffness_matrix_local2 ( long  eid,
long  ri,
long  ci,
matrix sm 
)
void transf_matrix ( ivector nodes,
matrix tmat 
)

The function assembles transformation matrix from local nodal coordinate system to the global coordinate system x_g = T x_l.

Parameters:
nodes - array containing node numbers
tmat - transformation matrix
Returns:
The function returns resulting transformation matrix in the parameter tmat.

PF, 20.12.2002, revised by JK, 28.11.2006

Definition at line 253 of file beamel3d.cpp.

References node::e1, node::e2, node::e3, fillm(), matrix::m, Mt, ivector::n, mechtop::nodes, and node::transf.

Referenced by initstr_matrix(), internal_forces(), load_matrix(), mass_matrix(), nodal_displ(), nodal_forces(), and stiffness_matrix().


Member Data Documentation

long* cncomp

array containing cumulative numbers of components of stress and strain tensors

Definition at line 56 of file beamel3d.h.

Referenced by beamel3d(), and ~beamel3d().

long intordism

order of integration of initial stress matrix

Definition at line 64 of file beamel3d.h.

Referenced by beamel3d().

long intordmm

order of integration of mass matrix

Definition at line 62 of file beamel3d.h.

Referenced by beamel3d(), load_matrix(), and mass_matrix().

long** intordsm

order of integration of stiffness matrix

Definition at line 60 of file beamel3d.h.

Referenced by beamel3d(), mechtop::give_intordsm(), and ~beamel3d().

long napfun

number of approximated functions on the element

Definition at line 58 of file beamel3d.h.

Referenced by beamel3d(), mechtop::give_napfun(), and load_matrix().

long nb

number of blocks

Definition at line 68 of file beamel3d.h.

Referenced by beamel3d(), mechtop::give_nb(), mechtop::give_nb_te(), and ~beamel3d().

long* ncomp

array containing numbers of components of stress and strain tensors

Definition at line 54 of file beamel3d.h.

Referenced by beamel3d(), and ~beamel3d().

long ndofe
long** nip

array of numbers of integration points in sets

Definition at line 66 of file beamel3d.h.

Referenced by beamel3d(), mechtop::give_nip(), and ~beamel3d().

long nne

stress/strain state

Definition at line 70 of file beamel3d.h.

Referenced by beamel3d(), and mechtop::give_ssst().

long tncomp

total number of components of stress and strain tensors

Definition at line 50 of file beamel3d.h.

Referenced by beamel3d(), mechtop::give_ncomp(), mechtop::give_tncomp(), initstr_matrix(), internal_forces(), load_matrix(), stiffness_matrix(), stiffness_matrix_local(), and stiffness_matrix_local2().

long tnip

total number of integration points

Definition at line 52 of file beamel3d.h.

Referenced by beamel3d(), and mechtop::give_tnip().


The documentation for this class was generated from the following files:

Generated by  doxygen 1.6.2