beamgen3d Class Reference

#include <beamgen3d.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)
 beamgen3d (void)
void bf_matrix (matrix &n, double s, double dl, double gy, double gz)
void ck_matrix (matrix &ck, double s, double c, double eh, double dl)
void geom_matrix (matrix &n, double s, double dl, double gy, double gz)
void internal_forces (long lcid, long eid, long ri, long ci, vector &ifor)
void internal_forces1 (long lcid, long eid, long ri, long ci, vector &ifor)
void load_matrix (long eid, matrix &lm)
void nodal_displ (long eid, long lcid)
void nodal_forces (long eid, long lcid)
void res_internal_forces (long lcid, long eid, vector &ifor)
void res_stiffness_matrix (long eid, matrix &sm)
void stiffness_matrix (long eid, long ri, long ci, matrix &sm)
void stiffness_matrixtor (matrix &sm, double dl, double e, double g, double gy, double gz, double *ixyz, double *ioyz)
void stiffness_matrixtor1 (matrix &sm, double dl, double e, double g, double gy, double gz, double *ixyz, double *ioyz)
void stiffness_matrixtor2 (matrix &sm, double dl, double e, double g, double gy, double gz, double *ixyz, double *ioyz, double *iro)
void transf_matrix (ivector &nodes, matrix &tmat)
 ~beamgen3d (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 beamgen3d describes threedimensional beam with general cross section where centre of shear may be different from center of gravity

PF, 2.10.2006

Definition at line 17 of file beamgen3d.h.


Constructor & Destructor Documentation

beamgen3d ( void   ) 

Definition at line 11 of file beamgen3d.cpp.

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

~beamgen3d ( void   ) 

Definition at line 44 of file beamgen3d.cpp.

References nb, and nip.


Member Function Documentation

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

function assembles transformation matrix from local to global system x_g = T x_l

columns of the transformation matrix are created by coordinates of local base vectors expressed in global 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 105 of file beamgen3d.cpp.

References fillm().

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

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.9.2006

Definition at line 311 of file beamgen3d.cpp.

References ll.

Referenced by load_matrix().

void ck_matrix ( matrix ck,
double  s,
double  c,
double  eh,
double  dl 
)

Definition at line 205 of file beamgen3d.cpp.

References fillm().

Referenced by internal_forces(), and stiffness_matrixtor().

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.9.2006

Definition at line 259 of file beamgen3d.cpp.

References fillm(), and ll.

void internal_forces ( long  lcid,
long  eid,
long  ri,
long  ci,
vector ifor 
)
void internal_forces1 ( long  lcid,
long  eid,
long  ri,
long  ci,
vector ifor 
)
void load_matrix ( long  eid,
matrix lm 
)
void nodal_displ ( long  eid,
long  lcid 
)

function stores nodal displacements and rotations in local 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

PF, 20.9.2006

Definition at line 877 of file beamgen3d.cpp.

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

Referenced by compute_ipstrains().

void nodal_forces ( long  eid,
long  lcid 
)

function computes nodal forces and moments expressed in local coodinate system this is equivalent to functions computing stresses 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 937 of file beamgen3d.cpp.

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

Referenced by compute_ipstresses(), and compute_nodestresses().

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

function computes internal forces

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

PF 20.9.2006

Definition at line 1012 of file beamgen3d.cpp.

References internal_forces1().

void res_stiffness_matrix ( long  eid,
matrix sm 
)

function computes stiffness matrix of 3D beam element

Parameters:
eid - number of element
sm - stiffness matrix

PF 20.9.2006

Definition at line 360 of file beamgen3d.cpp.

References matrix::m, Mp, stiffness_matrix(), and probdesc::zero.

Referenced by stiffmat().

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

function computes stiffness matrix of 3D beam element

Parameters:
eid - number of element
ri,ci - row and column indices
sm - stiffness matrix

PF, 20.9.2006

!! Mc->give_iomega (eid,ioyz);

Definition at line 381 of file beamgen3d.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(), glmatrixtransf(), element::ipp, ll, mechtop::locsystems(), mechmat::matstiff(), Mc, Mm, Mp, Mt, ndofe, nne, nodes, Out, stiffness_matrixtor(), tncomp, transf_matrix(), and probdesc::zero.

Referenced by nodal_forces(), and res_stiffness_matrix().

void stiffness_matrixtor ( matrix sm,
double  dl,
double  e,
double  g,
double  gy,
double  gz,
double *  ixyz,
double *  ioyz 
)

function computes consistent mass matrix of 3D beam element influence of inertial forces from rotations is accounted

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

PF, 20.9.2006

Definition at line 498 of file beamgen3d.cpp.

References ck_matrix().

Referenced by stiffness_matrix().

void stiffness_matrixtor1 ( matrix sm,
double  dl,
double  e,
double  g,
double  gy,
double  gz,
double *  ixyz,
double *  ioyz 
)

Definition at line 602 of file beamgen3d.cpp.

References ll.

void stiffness_matrixtor2 ( matrix sm,
double  dl,
double  e,
double  g,
double  gy,
double  gz,
double *  ixyz,
double *  ioyz,
double *  iro 
)

Definition at line 711 of file beamgen3d.cpp.

void transf_matrix ( ivector nodes,
matrix tmat 
)

function assembles transformation matrix from global to nodal system x_n = T x_g

Parameters:
nodes - array containing node numbers
tmat - transformation matrix

PF, 20.12.2002

Definition at line 63 of file beamgen3d.cpp.

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

Referenced by internal_forces(), internal_forces1(), load_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 53 of file beamgen3d.h.

Referenced by beamgen3d().

long intordism

order of integration of initial stress matrix

Definition at line 61 of file beamgen3d.h.

Referenced by beamgen3d().

long intordmm

order of integration of mass matrix

Definition at line 59 of file beamgen3d.h.

Referenced by beamgen3d(), and load_matrix().

long** intordsm

order of integration of stiffness matrix

Definition at line 57 of file beamgen3d.h.

Referenced by beamgen3d(), and mechtop::give_intordsm().

long napfun

number of approximated functions on the element

Definition at line 55 of file beamgen3d.h.

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

long nb

number of blocks

Definition at line 65 of file beamgen3d.h.

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

long* ncomp

array containing numbers of components of stress and strain tensors

Definition at line 51 of file beamgen3d.h.

Referenced by beamgen3d().

long ndofe

number of DOFs on the element

Definition at line 43 of file beamgen3d.h.

Referenced by beamgen3d(), mechtop::give_ndofe(), internal_forces(), internal_forces1(), load_matrix(), nodal_displ(), nodal_forces(), and stiffness_matrix().

long** nip

array of numbers of integration points in sets

Definition at line 63 of file beamgen3d.h.

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

long nne

number of nodes on one element

Definition at line 45 of file beamgen3d.h.

Referenced by beamgen3d(), mechtop::give_nne(), internal_forces(), internal_forces1(), load_matrix(), nodal_displ(), nodal_forces(), and stiffness_matrix().

stress/strain state

Definition at line 67 of file beamgen3d.h.

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

long tncomp

total number of components of stress and strain tensors

Definition at line 47 of file beamgen3d.h.

Referenced by beamgen3d(), mechtop::give_ncomp(), mechtop::give_tncomp(), internal_forces(), internal_forces1(), load_matrix(), and stiffness_matrix().

long tnip

total number of integration points

Definition at line 49 of file beamgen3d.h.

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


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

Generated by  doxygen 1.6.2