lintet Class Reference

#include <lintet.h>

List of all members.

Public Member Functions

double approx (vector &volcoord, vector &nodval)
double approx_nat (double xi, double eta, double zeta, vector &nodval)
void bd_matrix (long eid, long ri, long ci, matrix &bd)
void bf_matrix (matrix &n, double xi, double eta, double zeta)
void bf_matrix (matrix &n, vector &volcoord)
void compute_eigstress (long lcid, long eid, long ri, long ci)
double compute_error (long eid, double &e2, double &u2, double &sizel, vector *rsigfull)
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 dd_matrix (long eid, long ri, long ci, matrix &dd)
void eigstrain_forces (long eid, vector &nfor)
void eigstrain_forces (long lcid, long eid, long ri, long ci, vector &nfor, vector &x, vector &y, vector &z)
void elchar (long eid, matrix &spsig)
void elem_integration (integratedquant iq, long lcid, long eid, long ri, long ci, vector &nv, vector &x, vector &y, vector &z)
void elem_integration_quant (long eid, integratedquant iq, long lcid, vector &nv)
void geom_matrix (matrix &gm, vector &x, vector &y, vector &z, double &jac)
void geom_matrix_old (matrix &gm, vector &x, vector &y, vector &z)
void incr_internal_forces (long lcid, long eid, long ri, long ci, vector &ifor, vector &x, vector &y, vector &z)
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, vector &x, vector &y, vector &z)
void intpointval (long eid, vector &nodval, vector &ipval)
void ip_strains (long lcid, long eid, long ri, long ci)
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)
 lintet (void)
void load_matrix (long eid, matrix &lm)
void local_values (long lcid, long eid, long ri, long ci)
void locglob_nodeval (long is, vector &nv, double *tnv, vector &x, vector &y, vector &z)
void mass_matrix (long eid, matrix &mm)
void nod_eqother_ip (long eid)
void nod_strains_ip (long lcid, long eid, long ri, long ci)
void nod_stresses (long lcid, long eid, long ri, long ci)
void nod_stresses_ip (long lcid, long eid, long ri, long ci)
void node_forces_surf (long lcid, long eid, long *is, double *nv, vector &nf)
void node_forces_surf_old (long lcid, long eid, long *is, double *nv, vector &nf)
void nonloc_internal_forces (long lcid, long eid, long ri, long ci, vector &ifor, vector &x, vector &y, vector &z)
void ntdbr_vector (long eid, vector &ntdbr)
void ntn_matrix (long eid, matrix &ntn)
void res_eigstrain_forces (long eid, vector &nfor)
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_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 transf_matrix (ivector &nodes, matrix &tmat)
double volumeip (long eid, long ri, long ci)
 ~lintet (void)

Public Attributes

long * cncomp
 cumulative array of numbers of components of blocks
long intordb
 order of numerical interation on surfaces
long intordmm
 order of integration of mass matrix
long ** intordsm
 array of orders of integration of stiffness matrix
long napfun
 number of approximated functions on the element
long nb
 number of blocks
long * ncomp
 array of numbers of components of blocks
long ndofe
 number of DOFs on the element
long ned
 number of edges on one element
long ** nip
 array of numbers of integration points in sets
long nne
 number of nodes on one element
long nned
 number of nodes on one edge
long nnsurf
 number of nodes on one surface
long nsurf
 number of surfaces
strastrestate ssst
 stress/strain state
long tncomp
 total number of components of the strain and stress tensors
long tnip
 total number of integration points on element

Detailed Description

class lintet defines tetrahedral elements with linear approximation functions

JK

Definition at line 15 of file lintet.h.


Constructor & Destructor Documentation

lintet ( void   ) 

Definition at line 17 of file lintet.cpp.

References cncomp, intordb, intordmm, intordsm, napfun, nb, ncomp, ndofe, ned, nip, nne, nned, nnsurf, nsurf, spacestress, ssst, tncomp, and tnip.

~lintet ( void   ) 

Definition at line 78 of file lintet.cpp.

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


Member Function Documentation

double approx ( vector volcoord,
vector nodval 
)

function approximates function defined by nodal values

Parameters:
volcoord - volume coordinates
nodval - nodal values

20.8.2001

Definition at line 102 of file lintet.cpp.

References f, and scprd().

Referenced by interpolate(), intpointval(), ipcoord(), and mass_matrix().

double approx_nat ( double  xi,
double  eta,
double  zeta,
vector nodval 
)

function approximates function defined by nodal values

Parameters:
volcoord - volume coordinates
nodval - nodal values

20.8.2001

Definition at line 119 of file lintet.cpp.

References f, and scprd().

Referenced by inicipval(), and quadtet::intpointval2().

void bd_matrix ( long  eid,
long  ri,
long  ci,
matrix bd 
)

function computes B^T D matrix of one element this matrix is used in homogenization methods with prescribed stresses

Parameters:
eid - number of element
ri,ci - row and column indices
bd - B^T D matrix

17. 2. 2014

Definition at line 316 of file lintet.cpp.

References cmulm(), mechtop::elements, fillm(), geom_matrix(), mechtop::give_node_coord3d(), element::ipp, mechmat::matstiff(), Mm, Mt, mtxm(), ndofe, nne, and tncomp.

Referenced by stiffness_matrix().

void bf_matrix ( matrix n,
double  xi,
double  eta,
double  zeta 
)

Definition at line 150 of file lintet.cpp.

References fillm().

void bf_matrix ( matrix n,
vector volcoord 
)

function assembles matrix of base functions

Parameters:
n - matrix of base functions
volcoord - volume coordinates

24.8.2001

Definition at line 142 of file lintet.cpp.

References fillm().

Referenced by load_matrix(), mass_matrix(), node_forces_surf(), and node_forces_surf_old().

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

function computes correct eigen stresses caused by temperature at integration points on element

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

JK, 27.11.2006

Definition at line 1101 of file lintet.cpp.

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

Referenced by eigstrain_forces().

double compute_error ( long  eid,
double &  e2,
double &  u2,
double &  sizel,
vector rsigfull 
)

Function 1. 1. computes "e2" - square of energy norm of error of solution on element e2 = integral{ (sig_star - sig_roof)T * D_inv * (sig_star - sig_roof) } sig_star = recovered stress, values in nodes are defined by "rsigfull" and over element are interpolated by base functions sig_roof = stress obtained by FEM D_inv = inverse stiffness matrix of material 2. computes "u2" - square of energy norm of strain on element u2 = integral{ epsT * D * eps } eps = strain obtained by FEM D = stiffness matrix of material 3. computes area of element and adds to "volume" (sum of areas of previous elements) 4. computes "sizel" (characteristic size of element)

Parameters:
eid - element id
volume - sum of areas of previous elements
e2 - empty(returned) value;
u2 - empty(returned) value;
sizel - empty(returned) value;
rsigfull - 1Darray of stress in nodes; dimmension is ncomp*gt->nn after this manner {(val[0]; ... ;val[gt->nn])[0] ; ... ; (val[0]; ... ;val[gt->nn])[ncomp]}

created 3.5.2002, Ladislav Svoboda, termit@cml.fsv.cvut.cz

Definition at line 2098 of file lintet.cpp.

References vector::a, det3d(), mechtop::elements, fillv(), gauss_points_tet(), give_der_star(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), mechmat::givestrain(), invm(), element::ipp, mechmat::matstiff(), Mm, Mt, mxv(), mechtop::nn, nne, nodes, scprd(), subv(), tncomp, and vxmxv().

Referenced by adaptivity::compute_error().

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

TKo 7.2008

Definition at line 1012 of file lintet.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 increments of stresses at integration points on element

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

TKo 7.2008

Definition at line 1034 of file lintet.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

TKo, 7.2008

Definition at line 1079 of file lintet.cpp.

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

Referenced by nonloc_internal_forces().

void dd_matrix ( long  eid,
long  ri,
long  ci,
matrix dd 
)

function computes integral of D matrix of one element this matrix is used in homogenization methods with prescribed stresses

Parameters:
eid - number of element
ri,ci - row and column indices
bd - integral of D matrix

17. 2. 2014

Definition at line 351 of file lintet.cpp.

References vector::a, cmulm(), det3d(), mechtop::elements, fillm(), mechtop::give_node_coord3d(), element::ipp, mechmat::matstiff(), Mm, Mt, nne, and tncomp.

Referenced by stiffness_matrix().

void eigstrain_forces ( long  eid,
vector nfor 
)

function computes contributions from eigenstrains to the right hand side

Parameters:
eid - element id
ifor - vector of internal forces

JK, 22.4.2005

Definition at line 1968 of file lintet.cpp.

References cmulv(), cncomp, det3d(), mechtop::elements, fillv(), geom_matrix(), mechtop::give_node_coord3d(), mechmat::giveeigstrain(), element::ipp, mechmat::matstiff(), Mm, Mt, mtxv(), mxv(), vector::n, ncomp, ndofe, nne, and tncomp.

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

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,z - nodal coordinates

7.2008, TKo

Definition at line 864 of file lintet.cpp.

References compute_eigstress(), probdesc::eigstrcomp, eigstress, elem_integration(), and Mp.

Referenced by res_eigstrain_forces().

void elchar ( long  eid,
matrix spsig 
)

This function prepare array 'spder' for using in function `least_squarespr_default`. It allocates and fills up array 'spder' by derivatives(strain ro stress) in sampling points. For planeelemlt sampling points == main integration point == Gauss's point at the centroid of element.

Parameters:
eid - element id
spder - allocated array
flags - determines by which derivate is spder filled

created 3.5.2002, Ladislav Svoboda, termit@cml.fsv.cvut.cz

Definition at line 2167 of file lintet.cpp.

References vector::a, matrix::a, allocm(), mechtop::elements, mechmat::givestrain(), element::ipp, matrix::m, mechmat::matstiff(), Mm, Mt, mxv(), matrix::n, and tncomp.

Referenced by adaptivity::spr().

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

The function integrates selected stress type quantity over the finite element, i.e. it performs {} {B}^T {} d. It results in nodal values.

Parameters:
iq - type of integrated quantity (see alias.h)
lcid - number of load case
eid - element id
ri,ci - row and column indices
nv - nodal values (output)
x,y,z - node coordinates
Returns:
The function returns nodal values calculated in the vector nv

TKo 7.2008

Definition at line 1132 of file lintet.cpp.

References vector::a, addv(), cmulv(), cncomp, det3d(), mechtop::elements, fillv(), geom_matrix(), mechmat::givequantity(), element::ipp, Mm, Mt, mtxv(), ndofe, and tncomp.

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

void elem_integration_quant ( long  eid,
integratedquant  iq,
long  lcid,
vector nv 
)

The function integrates arbitrary selected quantity over the finite element, i.e. it performs {} {} d. It results in nodal values that can be used in the homogenization problems.

Parameters:
eid - element id
iq - type of integrated quantity (see alias.h)
lcid - number of load case
nv - nodal values (output)
Returns:
The function returns nodal values calculated in the vector nv

Created by TKo, 23.1.2015

Definition at line 1177 of file lintet.cpp.

References vector::a, cmulv(), cncomp, det3d(), mechtop::elements, fillv(), mechtop::give_node_coord3d(), mechmat::givequantity(), element::ipp, Mm, Mt, nne, and tncomp.

Referenced by elem_integration_quant().

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

function assembles geometric matrix vector of strains has following ordering eps=(e_xx, e_yy, e_zz, e_yz, e_zx, e_xy)

Parameters:
gm - geometric matrix
x,y,z - vectors of node coordinates

24.8.2001

Definition at line 208 of file lintet.cpp.

References vector::a, derivatives_3d(), dx_bf_lin_tet(), dy_bf_lin_tet(), and dz_bf_lin_tet().

Referenced by bd_matrix(), eigstrain_forces(), elem_integration(), ip_strains(), and stiffness_matrix().

void geom_matrix_old ( matrix gm,
vector x,
vector y,
vector z 
)

function assembles geometric matrix vector of strains has following ordering eps=(e_xx, e_yy, e_zz, e_yz, e_zx, e_xy)

Parameters:
gm - geometric matrix
x,y,z - vectors of node coordinates

24.8.2001

Definition at line 169 of file lintet.cpp.

References vector::a, det3d(), volb_3d(), volc_3d(), and vold_3d().

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

function computes increments of internal forces

Parameters:
lcid - load case id
eid - element id
ri,ci - row and column indices
ifor - vector of internal forces
x,y,z - vectors of nodal coordinates

TKo 7.2008

Definition at line 841 of file lintet.cpp.

References compute_nlstressincr(), elem_integration(), 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 1259 of file lintet.cpp.

References vector::a, allocv(), approx_nat(), destrv(), mechtop::elements, intpoints::eqother, gauss_points_tet(), 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,
vector x,
vector y,
vector z 
)

function computes internal forces

Parameters:
lcid - load case id
eid - element id
ri,ci - row and column indices
ifor - vector of internal forces
x,y,z - vectors of nodal coordinates

TKo 7.2008

Definition at line 793 of file lintet.cpp.

References compute_nlstress(), elem_integration(), and locstress.

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

Parameters:
eid - element id

JK, 22.4.2005

Definition at line 1934 of file lintet.cpp.

References approx().

Referenced by elem_intpointval().

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

function computes strains in main integration points of element

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

10.5.2002

Definition at line 506 of file lintet.cpp.

References allocm(), allocv(), copyv(), destrm(), destrv(), eldispl(), mechtop::elements, geom_matrix(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), element::ipp, lgvectortransf(), mechtop::locsystems(), Mm, Mt, mxv(), ndofe, nne, nodes, mechmat::storestrain(), tncomp, and transf_matrix().

Referenced by res_ip_strains().

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

function computes stresses in main integration points of element

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

10.5.2002

Definition at line 656 of file lintet.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 returns coordinates of integration points

Parameters:
eid - element id
ipp - integration point pointer
ri - row index
ci - column index
coord - vector of coordinates

19.1.2002

Definition at line 1211 of file lintet.cpp.

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

Referenced by ipcoord().

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

function computes volume appropriate to integration point

2.3.2004, JK

Definition at line 1373 of file lintet.cpp.

References vector::a, det3d(), mechtop::elements, mechtop::give_node_coord3d(), element::ipp, Mm, Mt, nne, and mechmat::storeipvol().

Referenced by ipvolume().

void load_matrix ( long  eid,
matrix lm 
)

function computes load matrix

Parameters:
eid - number of element
lm - load matrix

26.8.2001

Definition at line 447 of file lintet.cpp.

References matrix::a, bf_matrix(), det3d(), fillm(), gauss_points_tet(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), intordmm, matrix::m, Mt, matrix::n, napfun, ndofe, nne, nnj(), and nodes.

Referenced by loadmat().

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 1057 of file lintet.cpp.

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

Referenced by elem_local_values().

void locglob_nodeval ( long  is,
vector nv,
double *  tnv,
vector x,
vector y,
vector z 
)

Definition at line 1848 of file lintet.cpp.

References vector::a, matrix::a, cmulv(), mxv(), and scprd().

Referenced by node_forces_surf().

void mass_matrix ( long  eid,
matrix mm 
)

function computes mass matrix

Parameters:
eid - number of element
mm - mass matrix

26.8.2001

Definition at line 404 of file lintet.cpp.

References matrix::a, approx(), bf_matrix(), det3d(), fillm(), gauss_points_tet(), mechcrsec::give_density(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), intordmm, matrix::m, Mc, Mt, matrix::n, napfun, ndofe, nne, nnj(), and nodes.

Referenced by massmat().

void nod_eqother_ip ( long  eid  ) 
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

10.5.2002

Definition at line 547 of file lintet.cpp.

References cmulv(), mechtop::elements, mechtop::give_elemnodes(), mechmat::givestrain(), element::ipp, Mm, Mp, Mt, nne, nod, mechtop::nodes, node::storestrain(), probdesc::strainaver, tncomp, and volumeip().

Referenced by compute_nodestrains(), and strains().

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

function computes stresses in nodes

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

10.5.2002

Definition at line 677 of file lintet.cpp.

References cmulv(), mechtop::elements, mechtop::give_elemnodes(), mechmat::givestress(), element::ipp, Mm, Mp, Mt, nne, nod, mechtop::nodes, node::storestress(), probdesc::stressaver, tncomp, and volumeip().

Referenced by compute_nodestresses(), and stresses().

void node_forces_surf ( long  lcid,
long  eid,
long *  is,
double *  nv,
vector nf 
)

function computes nodal forces caused by presure on surface

Parameters:
lcid - number of load case
eid - element id
is - identification of surfaces
nv - nodal values
nf - nodal forces

12.12.2005, JK

Definition at line 1579 of file lintet.cpp.

References vector::a, addv(), allocm(), bf_matrix(), fillm(), fillv(), gauss_points_tr(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), intordb, jac2d_3d(), lgvectortransfblock(), locglob_nodeval(), mechtop::locsystems(), Mt, mxv(), napfun, ndofe, nne, nnj(), nodes, and transf_matrix().

Referenced by loadel::surfaceload().

void node_forces_surf_old ( long  lcid,
long  eid,
long *  is,
double *  nv,
vector nf 
)

function computes nodal forces caused by presure on surface

Parameters:
lcid - number of load case
eid - element id
is - identification of surfaces
nv - nodal values
nf - nodal forces

5.4.2005, JK pracuje pouze v globalnim souradnem systemu

Definition at line 1453 of file lintet.cpp.

References matrix::a, vector::a, addv(), bf_matrix(), gauss_points_tr(), mechtop::give_node_coord3d(), intordb, jac2d_3d(), matrix::m, Mt, mxv(), matrix::n, napfun, ndofe, nne, and nnj().

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

function computes internal forces for nonlocal models

Parameters:
lcid - load case id
eid - element id
ri,ci - row and column indices
ifor - vector of internal forces
x,y,z - vectors of nodal coordinates

TKo 7.2008

Definition at line 817 of file lintet.cpp.

References compute_nonloc_nlstress(), elem_integration(), and nonlocstress.

Referenced by res_nonloc_internal_forces().

void ntdbr_vector ( long  eid,
vector ntdbr 
)

Function integrates function "NT*D*B*r" (= NT*Stress) over whole element. N is matrix of basic functions. !!! Values in 'ntdbr' are stored unusually. Not after this manner {(val[0]; ... ;val[tncomp])[0] ; ... ; (val[0]; ... ;val[tncomp])[nne]} but after this manner {(val[0]; ... ;val[nne])[0] ; ... ; (val[0]; ... ;val[nne])[ntncomp]}

Parameters:
eid - element id
ntdbr - empty(returned) array, dimension is tncomp*nne

created 3.5.2002, Ladislav Svoboda, termit@cml.fsv.cvut.cz

Definition at line 2020 of file lintet.cpp.

References det3d(), mechtop::elements, fillv(), mechtop::give_node_coord3d(), mechmat::givestrain(), element::ipp, mechmat::matstiff(), Mm, Mt, mxv(), nne, and tncomp.

Referenced by z2_smoothing::compute_ntdbr().

void ntn_matrix ( long  eid,
matrix ntn 
)

Function integrates function "NT*N" over whole element. N is matrix of basic functions. !!! Matrix N is composed of three same blocks, it is computed only one third.

Parameters:
eid - element id
ntn - empty(returned) 2Darray, dimension is nne x nne

created 3.5.2002, Ladislav Svoboda, termit@cml.fsv.cvut.cz

Definition at line 2058 of file lintet.cpp.

References vector::a, det3d(), fillm(), mechtop::give_node_coord3d(), Mt, and nne.

Referenced by z2_smoothing::compute_ainv(), and z2_smoothing::compute_ntn_sky().

void res_eigstrain_forces ( long  eid,
vector nfor 
)

function computes contributions from eigenstrains to the right hand side

Parameters:
eid - element id
ifor - vector of internal forces

JK, 22.4.2005

Definition at line 1955 of file lintet.cpp.

References eigstrain_forces().

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

function computes resulting contributions from eigenstrains to the right hand side

Parameters:
lcid - load case id
eid - element id
nfor - vector of internal forces

TKo, 7.2008

Definition at line 981 of file lintet.cpp.

References copyv(), eigstrain_forces(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), 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 increment of internal forces

Parameters:
lcid - load case id
eid - element id
ifor - vector of internal forces

TKo, 7.2008

Definition at line 950 of file lintet.cpp.

References copyv(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), 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

Parameters:
lcid - load case id
eid - element id
ifor - vector of internal forces

TKo, 7.2008

Definition at line 888 of file lintet.cpp.

References copyv(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), 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 strains at integration points

Parameters:
lcid - load case id
eid - element id

TKo - 1.2009

Definition at line 489 of file lintet.cpp.

References ip_strains().

Referenced by compute_ipstrains(), and strains().

void res_ip_stresses ( long  lcid,
long  eid 
)

Definition at line 641 of file lintet.cpp.

References ip_stresses().

Referenced by compute_ipstresses(), and stresses().

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

function computes resulting internal forces for nonlocal models

Parameters:
lcid - load case id
eid - element id
ifor - vector of internal forces

TKo, 7.2008

Definition at line 919 of file lintet.cpp.

References copyv(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), 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 assembles resulting stiffness matrix of the element

Parameters:
eid - element id
sm - stiffness matrix

JK, 9.5.2002

Definition at line 380 of file lintet.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 one element

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

19.7.2001

Definition at line 285 of file lintet.cpp.

References vector::a, bdbjac(), det3d(), mechtop::elements, fillm(), geom_matrix(), mechtop::give_node_coord3d(), element::ipp, mechmat::matstiff(), Mm, Mt, ndofe, nne, and tncomp.

Referenced by res_stiffness_matrix().

void strains ( long  lcid,
long  eid,
long  ri,
long  ci 
)
void stresses ( long  lcid,
long  eid,
long  ri,
long  ci 
)
void transf_matrix ( ivector nodes,
matrix tmat 
)

function assembles transformation matrix

Parameters:
nodes - nodes of element
tmat - transformation matrix

Definition at line 246 of file lintet.cpp.

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

Referenced by ip_strains(), node_forces_surf(), res_eigstrain_forces(), res_incr_internal_forces(), res_internal_forces(), res_nonloc_internal_forces(), and res_stiffness_matrix().

double volumeip ( long  eid,
long  ri,
long  ci 
)

function computes volume appropriate to integration point

2.3.2004, JK

Definition at line 1393 of file lintet.cpp.

References vector::a, det3d(), mechtop::elements, mechtop::give_node_coord3d(), element::ipp, Mt, and nne.

Referenced by nod_eqother_ip(), nod_strains_ip(), and nod_stresses_ip().


Member Data Documentation

long* cncomp

cumulative array of numbers of components of blocks

Definition at line 118 of file lintet.h.

Referenced by eigstrain_forces(), elem_integration(), elem_integration_quant(), lintet(), and ~lintet().

long intordb

order of numerical interation on surfaces

Definition at line 110 of file lintet.h.

Referenced by lintet(), node_forces_surf(), and node_forces_surf_old().

long intordmm

order of integration of mass matrix

Definition at line 108 of file lintet.h.

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

long** intordsm

array of orders of integration of stiffness matrix

Definition at line 106 of file lintet.h.

Referenced by export_gid_gauss_pt(), mechtop::give_intordsm(), inicipval(), ipcoord(), lintet(), and ~lintet().

long napfun

number of approximated functions on the element

Definition at line 96 of file lintet.h.

Referenced by mechtop::give_napfun(), lintet(), load_matrix(), mass_matrix(), node_forces_surf(), and node_forces_surf_old().

long nb
long* ncomp

array of numbers of components of blocks

Definition at line 116 of file lintet.h.

Referenced by adaptivity::check_consistency(), eigstrain_forces(), lintet(), and ~lintet().

long ndofe
long ned

number of edges on one element

Definition at line 98 of file lintet.h.

Referenced by mechtop::give_ned(), and lintet().

long** nip

array of numbers of integration points in sets

Definition at line 112 of file lintet.h.

Referenced by adaptivity::check_consistency(), mechtop::give_nip(), lintet(), and ~lintet().

long nne
long nned

number of nodes on one edge

Definition at line 100 of file lintet.h.

Referenced by mechtop::give_nned(), and lintet().

long nnsurf

number of nodes on one surface

Definition at line 104 of file lintet.h.

Referenced by mechtop::give_nnsurf(), and lintet().

long nsurf

number of surfaces

Definition at line 102 of file lintet.h.

Referenced by mechtop::give_nsurf(), and lintet().

stress/strain state

Definition at line 120 of file lintet.h.

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

long tncomp
long tnip

total number of integration points on element

Definition at line 94 of file lintet.h.

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


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

Generated by  doxygen 1.6.2