barel2d Class Reference

#include <barel2d.h>

List of all members.

Public Member Functions

double approx (double xi, vector &nodval)
 barel2d (void)
void compute_eigstress (long lcid, long eid, long ri, long ci)
void compute_nlstress (long lcid, long eid, long ri, long ci)
void compute_nlstressincr (long lcid, long eid, long ri, long ci)
void compute_nonloc_nlstress (long lcid, long eid, long ri, long ci)
void eigstrain_forces (long lcid, long eid, long ri, long ci, vector &nfor)
void elem_integration (integratedquant iq, long lcid, long eid, long ri, long ci, vector &nv, vector &x, vector &y)
void geom_matrix (matrix &gm, vector &x, double &jac)
void geom_matrix (matrix &gm, double s, double c, double l)
void giveloccoord (vector &x, vector &y, vector &lx)
void incr_internal_forces (long lcid, long eid, long ri, long ci, vector &ifor)
void inicipval (long eid, long ri, long ci, matrix &nodval, inictype *ictn)
void internal_forces (long lcid, long eid, long ri, long ci, vector &ifor)
void intpointval (long eid, vector &nodval, vector &ipval)
void ip_elast_stresses (long lcid, long eid, long ri, long ci)
void ip_strains (long lcid, long eid, long ri, long ci, vector &r)
void ip_stresses (long lcid, long eid, long ri, long ci)
void ipcoord (long eid, long ipp, long ri, long ci, vector &coord)
void ipvolume (long eid, long ri, long ci)
void load_matrix (long eid, matrix &lm, vector &x)
void local_values (long lcid, long eid, long ri, long ci)
void mass_matrix (long eid, matrix &mm)
void mechq_nodval (long eid, vector &nodval, nontransquant qt)
void mechq_nodval_comp (long eid, vector &nodval, long ncnv, long nq, nontransquant *qt)
void nod_eqother_ip (long eid, long ri, long ci)
void nod_strains_comp (long lcid, long eid)
void nod_strains_ip (long lcid, long eid, long ri, long ci)
void nod_stresses_ip (long lcid, long eid, long ri, long ci)
void nonloc_internal_forces (long lcid, long eid, long ri, long ci, vector &ifor)
void res_eigstrain_forces (long lcid, long eid, vector &nfor)
void res_incr_internal_forces (long lcid, long eid, vector &ifor)
void res_internal_forces (long lcid, long eid, vector &ifor)
void res_ip_strains (long lcid, long eid)
void res_ip_stresses (long lcid, long eid)
void res_load_matrix (long eid, matrix &lm)
void res_mass_matrix (long eid, matrix &mm)
void res_nonloc_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 strains (long lcid, long eid, long ri, long ci)
void stresses (long lcid, long eid, long ri, long ci)
void tran_mat (vector &x, matrix &tran, vector &gx, vector &gy)
void transf_matrix (ivector &nodes, matrix &tmat)
 ~barel2d (void)

Public Attributes

long * cncomp
 array containing cumulative numbers of components of stress and strain tensors
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 blocks
long nne
 number of nodes on one element
strastrestate ssst
 stress/strain state
long tncomp
 total number of components of the strain and stress tensors
long tnip
 total number of integration point

Detailed Description

class barel2d defines onedimensional bar element with linear approximation functions

JK

Definition at line 16 of file barel2d.h.


Constructor & Destructor Documentation

barel2d ( void   ) 

Definition at line 13 of file barel2d.cpp.

References bar, cncomp, intordsm, napfun, nb, ncomp, ndofe, nip, nne, ssst, tncomp, and tnip.

~barel2d ( void   ) 

Definition at line 61 of file barel2d.cpp.

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


Member Function Documentation

double approx ( double  xi,
vector nodval 
)

Function approximates function defined by nodal values.

Parameters:
xi - natural coordinate
nodval - nodal values
Return values:
Function returns approximated value for the given natural coordinate xi.

JK

Definition at line 137 of file barel2d.cpp.

References vector::a, bf_lin_1d(), f, nne, and scprd().

Referenced by inicipval(), intpointval(), barelq3d::intpointval2(), barelq2d::intpointval2(), and ipcoord().

void compute_eigstress ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes correct stresses caused by eigenstrains at integration points on element.

Parameters:
lcid - number of load case
eid - element id
ri,ci - row and column indices

JK, 30.11.2006

Definition at line 1512 of file barel2d.cpp.

References mechmat::eigenstresses(), mechtop::elements, mechmat::giveeigstrain(), intpoints::grmid(), mechmat::ip, element::ipp, mechmat::matstiff(), Mm, Mt, mxv(), mechmat::storeeigstress(), and tncomp.

Referenced by eigstrain_forces().

void compute_nlstress ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes correct stresses at integration points on element.

Parameters:
lcid - number of load case
eid - element id
ri,ci - row and column indices

JK, 30.11.2006

Definition at line 1418 of file barel2d.cpp.

References mechmat::computenlstresses(), mechtop::elements, element::ipp, Mm, Mp, Mt, and probdesc::strcomp.

Referenced by elem_nlstresses(), and internal_forces().

void compute_nlstressincr ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes correct stress increments at integration points on element.

Parameters:
lcid - number of load case
eid - element id
ri,ci - row and column indices

JK, 30.11.2006

Definition at line 1441 of file barel2d.cpp.

References mechmat::computenlstressesincr(), mechtop::elements, element::ipp, Mm, Mp, Mt, and probdesc::strcomp.

Referenced by incr_internal_forces().

void compute_nonloc_nlstress ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes nonlocal correct stresses at integration points on element.

Parameters:
lcid - number of load case
eid - element id
ri,ci - row and column indices

JK, 30.11.2006

Definition at line 1489 of file barel2d.cpp.

References mechmat::compnonloc_nlstresses(), mechtop::elements, element::ipp, Mm, Mp, Mt, and probdesc::strcomp.

Referenced by nonloc_internal_forces().

void eigstrain_forces ( long  lcid,
long  eid,
long  ri,
long  ci,
vector nfor 
)

Function computes nodal forces caused by temperature changes.

Parameters:
eid - element id
ri,ci - row and column indices
nfor - array containing nodal forces
x,y - nodal coordinates
Return values:
Function returns nodal values of forces caused by eigenstrains in vector nfor.

30.11.2002, JK

Definition at line 1012 of file barel2d.cpp.

References compute_eigstress(), probdesc::eigstrcomp, eigstress, elem_integration(), mechtop::give_node_coord2d(), Mp, Mt, and nne.

Referenced by res_eigstrain_forces().

void elem_integration ( integratedquant  iq,
long  lcid,
long  eid,
long  ri,
long  ci,
vector nv,
vector x,
vector y 
)

Function integrates selected quantity over the finite element

Parameters:
iq - type of integrated quantity (see alias.h)
lcid - number of load case (used for certain quantity types)
eid - element id
ri,ci - row and column indices
nv - nodal values
x,y - node coordinates
Return values:
Function returns integrated values at nodes in vector nv.

JK, 30.11.2006

Definition at line 1557 of file barel2d.cpp.

References addv(), cmulv(), mechtop::elements, geom_matrix(), mechcrsec::give_area(), mechmat::givequantity(), element::ipp, Mc, Mm, Mt, mtxv(), ndofe, nne, nodes, print_err(), and tncomp.

Referenced by eigstrain_forces(), incr_internal_forces(), internal_forces(), and nonloc_internal_forces().

void geom_matrix ( matrix gm,
vector x,
double &  jac 
)

Function assembles strain-displacement (geometric) matrix. Geometric matrix is in local element coordinate system.

Parameters:
gm - geometric matrix
x - node coordinates
jac - Jacobian
Return values:
Function returns Jacobian in parameter jac and geometric matrix in parameter gm.

JK, 10.8.2001

Definition at line 113 of file barel2d.cpp.

References derivatives_1d(), and dx_bf_lin_1d().

void geom_matrix ( matrix gm,
double  s,
double  c,
double  l 
)

function assembles strain-displacement (geometric) matrix geometric matrix is in global coordinate system

Parameters:
gm - geometric matrix
s,c - sin and cos
l - length of the bar
Return values:
Function returns geometric matrix in parameter gm.

JK, 10.8.2001

Definition at line 90 of file barel2d.cpp.

Referenced by elem_integration(), ip_strains(), barelc::lower_cap_coup_matrix(), barelc::lower_cond_coup_matrix(), barelc::mainip_strains(), nod_strains_comp(), stiffness_matrix(), barelc::upper_cap_coup_matrix(), barelc::upper_cond_coup_matrix(), barelc::upper_cond_coup_vector(), barelc::upper_internal_forces(), and barelc::volume_rhs_vector().

void giveloccoord ( vector x,
vector y,
vector lx 
)

Function returns coordinates of nodes at local element coordinate system.

Parameters:
x - vector of nodal x-coordinates in global coordinate system
y - vector of nodal y-coordinates in global coordinate system
lx - vector of nodal coordinates at local element coordinate system
Return values:
Function returns vector of nodal coordinates at parameter lx.

Definition at line 491 of file barel2d.cpp.

References sqr.

Referenced by res_mass_matrix().

void incr_internal_forces ( long  lcid,
long  eid,
long  ri,
long  ci,
vector ifor 
)

Function computes increment of internal forces.

Parameters:
lcid - load case id
eid - element id
ri,ci - row and column indices
ifor - vector of internal forces
Return values:
Function returns nodal values of increments of internal forces in vector ifor.

JK, 29.4.2008

Definition at line 982 of file barel2d.cpp.

References compute_nlstressincr(), elem_integration(), mechtop::give_node_coord2d(), Mt, nne, and stressincr.

Referenced by res_incr_internal_forces().

void inicipval ( long  eid,
long  ri,
long  ci,
matrix nodval,
inictype ictn 
)

The function computes initial values of the given quantities at each integration point of the element from the nodal values given by the parameter nodval. Initial condition types must be the same for all nodes of the element.

Parameters:
eid - element id
ri - block row index
ci - block column index
nodval - nodal values of particular initial conditions. nodval[i][j] represents value of j-th initial condition at i-th node of the given element.
ictn - array of types of initial condition for each node of element. The type of initial condition determines which values are being specified in the node. (ictn[i] & inistrain) returns nonzero if nodal values of initial strains are specified (ictn[i] & inistress) returns nonzero if nodal values of initial stresses are specified (ictn[i] & iniother) returns nonzero if nodal values of initial values of eqother array are specified (ictn[i] & inicond) returns nonzero if nodal values of other initial conditions are specified

The function does not return anything.

Created by Tomas Koudelka 2004 Revised by Tomas Koudelka 03.2012

Definition at line 1302 of file barel2d.cpp.

References vector::a, allocv(), approx(), destrv(), mechtop::elements, intpoints::eqother, gauss_points(), mechmat::ic, inicond, iniother, inistrain, inistress, intordsm, mechmat::ip, element::ipp, memset(), Mm, Mt, matrix::n, nb, intpoints::ncompeqother, intpoints::ncompstr, nne, print_err(), intpoints::strain, and intpoints::stress.

Referenced by mechbclc::inicipval().

void internal_forces ( long  lcid,
long  eid,
long  ri,
long  ci,
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
Return values:
Function returns nodal values of internal forces in vector ifor.

12.8.2001

Definition at line 924 of file barel2d.cpp.

References compute_nlstress(), elem_integration(), mechtop::give_node_coord2d(), locstress, Mt, and nne.

Referenced by res_internal_forces().

void intpointval ( long  eid,
vector nodval,
vector ipval 
)

Function interpolates the nodal values to the integration points on the element quadratic approximation functions are used

Parameters:
eid - element id
nodval - nodal values
ipval - value at integration points
Return values:
Function returns approximated values at integration points in the vector ipval.

Definition at line 1177 of file barel2d.cpp.

References vector::a, allocv(), approx(), destrv(), mechtop::elements, gauss_points(), intordsm, element::ipp, Mt, and nb.

Referenced by elem_intpointval(), and elem_intpointval2().

void ip_elast_stresses ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes stress at integration points.

Parameters:
lcid - load case id
eid - element id
ri - row index
ci - column index

10.5.2002

Definition at line 762 of file barel2d.cpp.

References mechtop::elements, mechmat::givestrain(), intordsm, element::ipp, mechmat::matstiff(), Mm, Mt, mxv(), mechmat::storestress(), and tncomp.

void ip_strains ( long  lcid,
long  eid,
long  ri,
long  ci,
vector r 
)

Function computes strains at integration point.

Parameters:
lcid - load case id
eid - element id
ri - row index
ci - column index
r - nodal displacements

10.5.2002

Definition at line 548 of file barel2d.cpp.

References mechtop::elements, geom_matrix(), mechtop::give_node_coord2d(), element::ipp, Mm, Mt, mxv(), ndofe, nne, mechmat::storestrain(), and tncomp.

Referenced by res_ip_strains().

void ip_stresses ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes stress at integration points.

Parameters:
lcid - load case id
eid - element id
ri - row index
ci - column index

10.5.2002

Definition at line 739 of file barel2d.cpp.

References mechmat::computenlstresses(), mechtop::elements, element::ipp, Mm, Mp, Mt, and probdesc::strcomp.

Referenced by res_ip_stresses().

void ipcoord ( long  eid,
long  ipp,
long  ri,
long  ci,
vector coord 
)

Function assembles coordinates of integration points.

Parameters:
eid - element id
ipp - integration point pointer
ri - row index
ci - column index
coord - array containing coordinates of integration points
Return values:
Function returns coordinates of integration points in the vector coord.

TKo, 7.2008

Definition at line 1222 of file barel2d.cpp.

References vector::a, approx(), mechtop::elements, gauss_points(), mechtop::give_node_coord2d(), intordsm, element::ipp, Mt, and nne.

void ipvolume ( long  eid,
long  ri,
long  ci 
)

Function computes volume appropriate to integration point.

Parameters:
eid - element id
ri - row index of the given integration point block
ci - column index of the given integration point block

Definition at line 1255 of file barel2d.cpp.

References mechtop::elements, mechcrsec::give_area(), mechtop::give_node_coord2d(), element::ipp, Mc, Mm, Mt, nne, print_err(), and mechmat::storeipvol().

void load_matrix ( long  eid,
matrix lm,
vector x 
)

Function computes load matrix.

Parameters:
eid - element id
lm - load matrix
x - vector of node coordinates
Return values:
Function returns load matrix in parameter lm.

JK, 12.10.2008

Definition at line 400 of file barel2d.cpp.

References fillm(), mechcrsec::give_area(), and Mc.

Referenced by res_load_matrix().

void local_values ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes local values which will be used for averageing in nonlocal models at integration points. Mp->nonlocphase have to be 1.

Parameters:
lcid - number of load case
eid - element id
ri,ci - row and column indices

TKo, 7.2008

Definition at line 1465 of file barel2d.cpp.

References mechmat::computenlstresses(), mechtop::elements, element::ipp, Mm, Mp, Mt, and probdesc::strcomp.

Referenced by elem_local_values().

void mass_matrix ( long  eid,
matrix mm 
)

Function computes mass matrix of twodimensional bar element.

Parameters:
eid - element id
mm - mass matrix
Return values:
Function returns mass matrix at parameter mm.

JK, 10.8.2001

Definition at line 312 of file barel2d.cpp.

References fillm(), mechcrsec::give_area(), mechcrsec::give_densitye(), mechtop::give_node_coord2d(), Mc, Mt, and nne.

Referenced by massmat(), and res_mass_matrix().

void mechq_nodval ( long  eid,
vector nodval,
nontransquant  qt 
)

Function approximates nonmechanical quantities in nodes of element. The nodal values are taken from the closest integration point.

Parameters:
eid - element id
nodval - vector of nodal values
qt - type of mechanical quantity
Returns:
The function does not return anything.

12/06/2012 TKr

Definition at line 1610 of file barel2d.cpp.

References mechtop::elements, mechmat::givemechq(), intordsm, element::ipp, Mm, Mt, nne, and nodip_bar().

Referenced by elem_mechq_nodval().

void mechq_nodval_comp ( long  eid,
vector nodval,
long  ncnv,
long  nq,
nontransquant qt 
)

Function computes mechanical quantities in nodes of element.

Parameters:
eid - element id
nodval - vector of nodal values of all required quantities, i.e., nodal value of i-th quantity in j-th node is given by nodval[i*ncnv+j] where ncnv is the number of calculated nodes on eid-th element.
ncnv - number of computed nodes on element (only first ncnv of nodes is calculated)
nq - number of required mechanical quantities
qt - array of types of required mechanical quantities
Returns:
The function does not return anything.

Created by Tomas Koudelka, 29.11.2013

Definition at line 1645 of file barel2d.cpp.

References mechmat::compnonloc_nlstresses(), mechmat::computenlstresses(), intpoints::copy(), mechtop::elements, intpoints::eqother, mechtop::give_elemnodes(), mechmat::give_num_averq(), mechmat::givemechq(), mechmat::givenonlocid(), intordsm, mechmat::ip, element::ipp, local, probdesc::matmodel, Mb, Mm, Mp, Mt, intpoints::ncompeqother, intpoints::ncompother, intpoints::ncompstr, mechbclc::nlc, nne, nod_strains_comp(), mechtop::nodes, nodip_bar(), intpoints::nonloc, nonlocal, intpoints::other, print_err(), mechmat::storeeqother(), mechmat::storenonloc(), mechmat::storeother(), mechmat::storestrain(), mechmat::storestress(), intpoints::strain, node::strain, intpoints::stress, and mechmat::updateipvalmat().

Referenced by elem_mechq_nodval_comp().

void nod_eqother_ip ( long  eid,
long  ri,
long  ci 
)

Function computes other values in nodes of element.

Parameters:
eid - element id

10.5.2002

Definition at line 882 of file barel2d.cpp.

References vector::a, allocv(), destrv(), mechtop::elements, mechtop::give_elemnodes(), mechmat::giveeqother(), mechmat::givencompeqother(), intordsm, element::ipp, Mm, Mt, nne, nod, mechtop::nodes, nodip_bar(), and node::storeother().

Referenced by compute_nodeothers().

void nod_strains_comp ( long  lcid,
long  eid 
)

The function computes nodal strains directly.

Parameters:
lcid - load case id
eid - element id

Tomas Koudelka, 14.11.2013

Definition at line 622 of file barel2d.cpp.

References vector::a, eldispl(), geom_matrix(), mechtop::give_elemnodes(), mechtop::give_node_coord2d(), Mt, mxv(), ndofe, nne, mechtop::nodes, node::storestrain(), and tncomp.

Referenced by mechq_nodval_comp(), and nodestrains_comp().

void nod_strains_ip ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes strains in nodes of element.

Parameters:
lcid - load case id
eid - element id
ri,ci - row and column indices

10.5.2002

Definition at line 588 of file barel2d.cpp.

References mechtop::elements, mechtop::give_elemnodes(), mechmat::givestrain(), intordsm, element::ipp, Mm, Mt, nne, nod, mechtop::nodes, nodip_bar(), node::storestrain(), and tncomp.

Referenced by compute_nodestrains(), and strains().

void nod_stresses_ip ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes stresses at nodes of element.

Parameters:
lcid - load case id
eid - element id
ri,ci - row and column indices

10.5.2002

Definition at line 795 of file barel2d.cpp.

References mechtop::elements, mechtop::give_elemnodes(), mechmat::givestress(), intordsm, element::ipp, Mm, Mt, nne, nod, mechtop::nodes, nodip_bar(), node::storestress(), and tncomp.

Referenced by compute_nodestresses(), and stresses().

void nonloc_internal_forces ( long  lcid,
long  eid,
long  ri,
long  ci,
vector ifor 
)

Function computes nonlocal internal forces.

Parameters:
lcid - number of load case
eid - element id
ri,ci - row and column indices
ifor - vector of internal forces
Return values:
Function returns nodal values of nonlocal internal forces in vector ifor.

TKo, 7.2008

Definition at line 955 of file barel2d.cpp.

References compute_nonloc_nlstress(), elem_integration(), mechtop::give_node_coord2d(), Mt, nne, and nonlocstress.

Referenced by res_nonloc_internal_forces().

void res_eigstrain_forces ( long  lcid,
long  eid,
vector nfor 
)

Function computes resulting internal forces caused by eigenstrains.

Parameters:
lcid - load case id
eid - element id
nfor - vector of nodal forces caused by eigenstrains
Return values:
Function returns resulting increments of internal forces in parameter nfor.

12.8.2001

Definition at line 1144 of file barel2d.cpp.

References copyv(), eigstrain_forces(), mechtop::give_elemnodes(), glvectortransf(), mechtop::locsystems(), Mt, ndofe, nne, nodes, and transf_matrix().

Referenced by elem_eigstrain_forces().

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

Function computes resulting increments of internal forces.

Parameters:
lcid - load case id
eid - element id
ifor - vector of internal forces at nodes
Return values:
Function returns resulting increments of internal forces in parameter ifor.

29.4.2008

Definition at line 1110 of file barel2d.cpp.

References copyv(), mechtop::give_elemnodes(), glvectortransf(), incr_internal_forces(), mechtop::locsystems(), Mt, ndofe, nne, nodes, and transf_matrix().

Referenced by elem_incr_internal_forces().

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

Function computes resulting internal forces. If it is required, nodal values are transformed to the local coordinate systems.

Parameters:
lcid - load case id
eid - element id
ifor - vector of internal forces at nodes
Return values:
Function returns resulting internal forces in parameter ifor.

12.8.2001

Definition at line 1042 of file barel2d.cpp.

References copyv(), mechtop::give_elemnodes(), glvectortransf(), internal_forces(), mechtop::locsystems(), Mt, ndofe, nne, nodes, and transf_matrix().

Referenced by elem_internal_forces().

void res_ip_strains ( long  lcid,
long  eid 
)

Function computes strain at integration points.

Parameters:
lcid - load case id
eid - element id

JK

Definition at line 511 of file barel2d.cpp.

References copyv(), eldispl(), mechtop::give_elemnodes(), ip_strains(), lgvectortransf(), mechtop::locsystems(), Mt, ndofe, nne, nodes, and transf_matrix().

Referenced by compute_ipstrains().

void res_ip_stresses ( long  lcid,
long  eid 
)

Function computes stress at integration points.

Parameters:
lcid - load case id
eid - element id

JK

Definition at line 722 of file barel2d.cpp.

References ip_stresses().

Referenced by compute_ipstresses().

void res_load_matrix ( long  eid,
matrix lm 
)

Function computes load matrix of the plane stress rectangular finite element with bilinear approximation functions. Load vector is obtained after premultiplying load matrix by nodal load values.

Parameters:
eid - number of element
lm - load matrix

JK, 25.7.2001

Definition at line 450 of file barel2d.cpp.

References mechtop::give_elemnodes(), mechtop::give_node_coord2d(), glmatrixtransf(), lgmatrixtransfblock(), load_matrix(), mechtop::locsystems(), Mt, ndofe, nne, nodes, tran_mat(), and transf_matrix().

Referenced by loadmat().

void res_mass_matrix ( long  eid,
matrix mm 
)

Function computes mass matrix of one element. If it is required, nodal values are transformed to the local coordinate systems.

Parameters:
eid - number of element
mm - mass matrix
Return values:
Function returns resulting mass matrix in parameter mm

Definition at line 365 of file barel2d.cpp.

References mechtop::give_elemnodes(), mechtop::give_node_coord2d(), giveloccoord(), glmatrixtransf(), mechtop::locsystems(), mass_matrix(), Mt, ndofe, nne, nodes, and transf_matrix().

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

Function computes resulting increment of nonlocal internal forces.

Parameters:
lcid - load case id
eid - element id
ifor - vector of internal forces at nodes
Return values:
Function returns resulting nonlocal internal forces in parameter ifor.

TKo 7.2008

Definition at line 1076 of file barel2d.cpp.

References copyv(), mechtop::give_elemnodes(), glvectortransf(), mechtop::locsystems(), Mt, ndofe, nne, nodes, nonloc_internal_forces(), and transf_matrix().

Referenced by elem_nonloc_internal_forces().

void res_stiffness_matrix ( long  eid,
matrix sm 
)

Function computes stiffness matrix of one element. If it is required, nodal values are transformed to the local coordinate systems.

Parameters:
eid - element id
sm - stiffness matrix
Return values:
Function returns resulting stiffness matrix in parameter sm

JK, 10.8.2001

Definition at line 282 of file barel2d.cpp.

References mechtop::give_elemnodes(), glmatrixtransf(), mechtop::locsystems(), Mt, ndofe, nne, nodes, stiffness_matrix(), and transf_matrix().

Referenced by stiffmat().

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

Function computes stiffness matrix of twodimensional bar element.

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

JK, 10.8.2001

Definition at line 234 of file barel2d.cpp.

References matrix::a, bdbj(), mechtop::elements, fillm(), geom_matrix(), mechcrsec::give_area(), mechtop::give_node_coord2d(), element::ipp, mechmat::matstiff(), Mc, Mm, Mt, ndofe, nne, nodes, print_err(), and tncomp.

Referenced by res_stiffness_matrix().

void strains ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes strains at strain points.

Parameters:
lcid - load case id
eid - element id
ri,ci - row and column indices

JK

Definition at line 668 of file barel2d.cpp.

References allocv(), destrv(), enodes, aepoints::give_aepcoord(), aepoints::give_naep(), aepoints::give_ncomp(), aepoints::give_sid(), intpts, Mm, Mp, nod_strains_ip(), nowhere, print_err(), aepoints::storevalues(), mechmat::stra, probdesc::strainaver, aepoints::tape, and userdefined.

void stresses ( long  lcid,
long  eid,
long  ri,
long  ci 
)

Function computes stresses at required positions.

Parameters:
lcid - load case id
eid - element id
ri,ci - row and column indices

Definition at line 830 of file barel2d.cpp.

References allocv(), destrv(), enodes, aepoints::give_aepcoord(), aepoints::give_naep(), aepoints::give_ncomp(), aepoints::give_sid(), intpts, Mm, Mp, nod_stresses_ip(), nowhere, print_err(), aepoints::storevalues(), mechmat::stre, probdesc::stressaver, aepoints::tape, and userdefined.

Referenced by computestresses().

void tran_mat ( vector x,
matrix tran,
vector gx,
vector gy 
)

Function assembles transformation matrix T between coordinate system defined on element and global coordinate system.

T x_local = x_global x_local = T^T x_global

direction vector denotes vector defining local coordinate axis

Parameters:
x - vector of nodal coordinates in element coordinate system
tran - tranformation matrix T
gx,gy - vectors of nodal coordinates in global coordinate systems

JK, 12.10.2008

Definition at line 194 of file barel2d.cpp.

References fillm().

Referenced by res_load_matrix().

void transf_matrix ( ivector nodes,
matrix tmat 
)

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

JK, 3.2.2002

Definition at line 160 of file barel2d.cpp.

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

Referenced by res_eigstrain_forces(), res_incr_internal_forces(), res_internal_forces(), res_ip_strains(), res_load_matrix(), res_mass_matrix(), res_nonloc_internal_forces(), and res_stiffness_matrix().


Member Data Documentation

long* cncomp

array containing cumulative numbers of components of stress and strain tensors

Definition at line 83 of file barel2d.h.

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

long intordmm

order of integration of mass matrix

Definition at line 89 of file barel2d.h.

long** intordsm
long napfun

number of approximated functions on the element

Definition at line 85 of file barel2d.h.

Referenced by barel2d(), and mechtop::give_napfun().

long nb

number of blocks

Definition at line 93 of file barel2d.h.

Referenced by barel2d(), barelc::barelc(), mechtop::give_nb(), mechtop::give_nb_te(), inicipval(), intpointval(), and ~barel2d().

long* ncomp

array containing numbers of components of stress and strain tensors

Definition at line 81 of file barel2d.h.

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

long ndofe
long** nip

array of numbers of integration points in blocks

Definition at line 91 of file barel2d.h.

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

long nne

stress/strain state

Definition at line 95 of file barel2d.h.

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

long tncomp
long tnip

total number of integration point

Definition at line 97 of file barel2d.h.

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


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

Generated by  doxygen 1.6.2