#include <lintetrot.h>
Public Member Functions | |
double | approx (vector &volcoord, vector &nodval) |
double | approx_nat (double xi, double eta, double zeta, vector &nodval) |
void | bf_matrix (matrix &n, double xi, double eta, double zeta) |
void | bf_matrix (matrix &n, vector &volcoord, vector &x, vector &y, vector &z) |
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 | dmatblock (matrix &dd, matrix &d) |
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 | elem_integration (integratedquant iq, long lcid, long eid, long ri, long ci, vector &nv, vector &x, vector &y, vector &z) |
void | geom_matrix (matrix &gm, vector &x, vector &y, vector &z, vector &volcoord) |
void | geom_matrix_shear (matrix &gm, vector &x, vector &y, vector &z, vector &volcoord) |
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, vector &x, vector &y, vector &z, vector &r) |
void | ipcoord (long eid, long ipp, long ri, long ci, vector &coord) |
void | ipvolume (long eid, long ri, long ci) |
lintetrot (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 lcid, 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 | 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 | tran_side (matrix &t12, matrix &t13, matrix &t14, matrix &t23, matrix &t34, matrix &t42, vector &x, vector &y, vector &z) |
void | transf_matrix (ivector &nodes, matrix &tmat) |
double | volumeip (long eid, long ri, long ci) |
~lintetrot (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 |
class lintetrot defines tetrahedral elements with linear approximation functions and rotational degrees of freedom
PF, 2.10.2008
Definition at line 16 of file lintetrot.h.
lintetrot | ( | void | ) |
~lintetrot | ( | void | ) |
function approximates function defined by nodal values
volcoord | - volume coordinates | |
nodval | - nodal values |
PF, 2.10.2008
Definition at line 100 of file lintetrot.cpp.
Referenced by intpointval(), ipcoord(), and mass_matrix().
double approx_nat | ( | double | xi, | |
double | eta, | |||
double | zeta, | |||
vector & | nodval | |||
) |
function approximates function defined by nodal values
volcoord | - volume coordinates | |
nodval | - nodal values |
20.8.2001
Definition at line 117 of file lintetrot.cpp.
Referenced by inicipval().
void bf_matrix | ( | matrix & | n, | |
double | xi, | |||
double | eta, | |||
double | zeta | |||
) |
Definition at line 259 of file lintetrot.cpp.
References fillm().
function assembles matrix of base functions
n | - matrix of base functions | |
volcoord | - volume coordinates |
2.10.2008
Definition at line 214 of file lintetrot.cpp.
References fillm(), and tran_side().
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
lcid | - number of load case | |
eid | - element id | |
ri,ci | - row and column indices |
JK, 27.11.2006
Definition at line 1243 of file lintetrot.cpp.
References mechtop::elements, mechmat::giveeigstrain(), 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
lcid | - number of load case | |
eid | - element id | |
ri,ci | - row and column indices |
TKo 7.2008
Definition at line 1152 of file lintetrot.cpp.
References mechmat::computenlstresses(), mechtop::elements, intordsm, element::ipp, Mm, Mp, Mt, and probdesc::strcomp.
Referenced by elem_nlstresses(), internal_forces(), and res_ip_stresses().
void compute_nlstressincr | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci | |||
) |
function computes correct increments of stresses at integration points on element
lcid | - number of load case | |
eid | - element id | |
ri,ci | - row and column indices |
TKo 7.2008
Definition at line 1176 of file lintetrot.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
lcid | - number of load case | |
eid | - element id | |
ri,ci | - row and column indices |
TKo, 7.2008
Definition at line 1221 of file lintetrot.cpp.
References mechmat::compnonloc_nlstresses(), mechtop::elements, element::ipp, Mm, Mp, Mt, and probdesc::strcomp.
Referenced by nonloc_internal_forces().
function assembles stiffness matrices for bending and shear
function extracts components of the bending part of stiffness matrix of the material to the matrix db
10.8.2008
Definition at line 505 of file lintetrot.cpp.
References intordsm.
Referenced by stiffness_matrix().
void eigstrain_forces | ( | long | eid, | |
vector & | nfor | |||
) |
function computes contributions from eigenstrains to the right hand side
eid | - element id | |
ifor | - vector of internal forces |
JK, 22.4.2005
Definition at line 2049 of file lintetrot.cpp.
References cmulv(), cncomp, det3d(), mechtop::elements, fillv(), 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
lcid | - load case id | |
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 1004 of file lintetrot.cpp.
References compute_eigstress(), probdesc::eigstrcomp, eigstress, elem_integration(), and Mp.
Referenced by res_eigstrain_forces().
void elem_integration | ( | integratedquant | iq, | |
long | lcid, | |||
long | eid, | |||
long | ri, | |||
long | ci, | |||
vector & | nv, | |||
vector & | x, | |||
vector & | y, | |||
vector & | z | |||
) |
function integrates selected quantity over the finite element it results in nodal values
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 | |
x,y,z | - node coordinates |
TKo 7.2008
Definition at line 1272 of file lintetrot.cpp.
References vector::a, addv(), cmulv(), cncomp, det3d(), mechtop::elements, fillv(), mechmat::givequantity(), element::ipp, Mm, Mt, mtxv(), ndofe, and tncomp.
Referenced by eigstrain_forces(), incr_internal_forces(), internal_forces(), and nonloc_internal_forces().
function assembles geometric (strain/displacement) matrix vector of strains has following ordering eps=(e_xx, e_yy, e_zz, e_yz, e_zx, e_xy)
gm | - geometric matrix | |
x,y,z | - vectors of node coordinates | |
volcoord | - array of volume coordinates |
PF, 2.10.2008
Definition at line 279 of file lintetrot.cpp.
References vector::a, det3d(), nne, tran_side(), volb_3d(), volc_3d(), and vold_3d().
Referenced by ip_strains(), and stiffness_matrix().
function assembles shear part of geometric (strain/displacement) matrix
gm | - geometric matrix | |
x,y,z | - vectors of nodal coordinates | |
volcoord | - array of volume coordinates |
PF, 2.10.2008
Definition at line 404 of file lintetrot.cpp.
References vector::a, det3d(), nne, tran_side(), volb_3d(), volc_3d(), and vold_3d().
Referenced by stiffness_matrix().
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
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 980 of file lintetrot.cpp.
References compute_nlstressincr(), elem_integration(), and stressincr.
Referenced by res_incr_internal_forces().
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.
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 1361 of file lintetrot.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.
void internal_forces | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci, | |||
vector & | ifor, | |||
vector & | x, | |||
vector & | y, | |||
vector & | z | |||
) |
function computes internal forces
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 932 of file lintetrot.cpp.
References compute_nlstress(), elem_integration(), and locstress.
Referenced by res_internal_forces().
function interpolates the nodal values to the integration points on the element
eid | - element id |
JK, 22.4.2005
Definition at line 2015 of file lintetrot.cpp.
References approx().
Referenced by elem_intpointval().
void ip_strains | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci, | |||
vector & | x, | |||
vector & | y, | |||
vector & | z, | |||
vector & | r | |||
) |
function computes strains in main integration points of element
lcid | - load case id | |
eid | - element id | |
ri | - row index | |
ci | - column index | |
x,y,z | - vectors of node coordinates | |
r | - vector of node displacements |
JK, 2.10.2008
Definition at line 784 of file lintetrot.cpp.
References vector::a, allocv(), destrv(), mechtop::elements, gauss_points_tet(), geom_matrix(), intordsm, element::ipp, Mm, Mt, mxv(), ndofe, mechmat::storestrain(), and tncomp.
Referenced by res_ip_strains().
void ipcoord | ( | long | eid, | |
long | ipp, | |||
long | ri, | |||
long | ci, | |||
vector & | coord | |||
) |
function returns coordinates of integration points
eid | - element id | |
ipp | - integration point pointer | |
ri | - row index | |
ci | - column index | |
coord | - vector of coordinates |
19.1.2002
Definition at line 1313 of file lintetrot.cpp.
References vector::a, approx(), mechtop::elements, gauss_points_tet(), mechtop::give_node_coord3d(), intordsm, element::ipp, Mt, and nne.
void ipvolume | ( | long | eid, | |
long | ri, | |||
long | ci | |||
) |
function computes volume appropriate to integration point
2.3.2004, JK
Definition at line 1477 of file lintetrot.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
eid | - number of element | |
lm | - load matrix |
26.8.2001
Definition at line 710 of file lintetrot.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.
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.
lcid | - number of load case | |
eid | - element id | |
ri,ci | - row and column indices |
TKo 7.2008
Definition at line 1199 of file lintetrot.cpp.
References mechmat::computenlstresses(), mechtop::elements, element::ipp, Mm, Mp, Mt, and probdesc::strcomp.
Referenced by elem_local_values().
Definition at line 1929 of file lintetrot.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
eid | - number of element | |
mm | - mass matrix |
26.8.2001
Definition at line 667 of file lintetrot.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 | lcid, | |
long | eid | |||
) |
Definition at line 1512 of file lintetrot.cpp.
References vector::a, allocv(), cmulv(), destrv(), mechtop::elements, mechtop::give_elemnodes(), mechmat::giveeqother(), mechmat::givencompeqother(), element::ipp, Mm, Mp, Mt, nne, nod, mechtop::nodes, probdesc::otheraver, node::storeother(), and volumeip().
void nod_strains_ip | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci | |||
) |
function computes strains in nodes of element
lcid | - load case id | |
eid | - element id |
10.5.2002
Definition at line 820 of file lintetrot.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().
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
lcid | - load case id | |
eid | - element id | |
ri,ci | - row and column indices |
10.5.2002
Definition at line 879 of file lintetrot.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().
void node_forces_surf | ( | long | lcid, | |
long | eid, | |||
long * | is, | |||
double * | nv, | |||
vector & | nf | |||
) |
function computes nodal forces caused by presure on surface
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 1681 of file lintetrot.cpp.
References matrix::a, vector::a, addv(), bf_matrix(), gauss_points_tr(), mechtop::give_node_coord3d(), intordb, jac2d_3d(), locglob_nodeval(), matrix::m, Mt, mxv(), matrix::n, napfun, ndofe, nne, and nnj().
void node_forces_surf_old | ( | long | lcid, | |
long | eid, | |||
long * | is, | |||
double * | nv, | |||
vector & | nf | |||
) |
function computes nodal forces caused by presure on surface
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 1555 of file lintetrot.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
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 956 of file lintetrot.cpp.
References compute_nonloc_nlstress(), elem_integration(), and nonlocstress.
Referenced by res_nonloc_internal_forces().
void res_eigstrain_forces | ( | long | eid, | |
vector & | nfor | |||
) |
function computes contributions from eigenstrains to the right hand side
eid | - element id | |
ifor | - vector of internal forces |
JK, 22.4.2005
Definition at line 2036 of file lintetrot.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
lcid | - load case id | |
eid | - element id | |
nfor | - vector of internal forces |
TKo, 7.2008
Definition at line 1121 of file lintetrot.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
lcid | - load case id | |
eid | - element id | |
ifor | - vector of internal forces |
TKo, 7.2008
Definition at line 1090 of file lintetrot.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
lcid | - load case id | |
eid | - element id | |
ifor | - vector of internal forces |
TKo, 7.2008
Definition at line 1028 of file lintetrot.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 | |||
) |
Definition at line 746 of file lintetrot.cpp.
References allocm(), allocv(), copyv(), destrm(), destrv(), eldispl(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), ip_strains(), lgvectortransf(), mechtop::locsystems(), Mt, ndofe, nne, nodes, tncomp, and transf_matrix().
Referenced by compute_ipstrains().
void res_ip_stresses | ( | long | lcid, | |
long | eid | |||
) |
function computes stresses at integration points
lcid | - load case id | |
eid | - element id |
JK, 2.10.2008
Definition at line 865 of file lintetrot.cpp.
References compute_nlstress().
Referenced by compute_ipstresses().
void res_nonloc_internal_forces | ( | long | lcid, | |
long | eid, | |||
vector & | ifor | |||
) |
function computes resulting internal forces for nonlocal models
lcid | - load case id | |
eid | - element id | |
ifor | - vector of internal forces |
TKo, 7.2008
Definition at line 1059 of file lintetrot.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
eid | - element id | |
sm | - stiffness matrix |
PF, 2.10.2008
Definition at line 642 of file lintetrot.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
eid | - number of element | |
ri,ci | - row and column indices | |
sm | - stiffness matrix |
PF, 2.10.2008
Definition at line 558 of file lintetrot.cpp.
References vector::a, allocm(), allocv(), bdbjac(), destrm(), destrv(), det3d(), dmatblock(), mechtop::elements, fillm(), gauss_points_tet(), geom_matrix(), geom_matrix_shear(), mechtop::give_node_coord3d(), intordsm, element::ipp, mechmat::matstiff(), Mm, Mt, ncomp, ndofe, nne, and tncomp.
Referenced by res_stiffness_matrix().
void strains | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci | |||
) |
Definition at line 851 of file lintetrot.cpp.
void stresses | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci | |||
) |
Definition at line 910 of file lintetrot.cpp.
Referenced by computestresses().
void tran_side | ( | matrix & | t12, | |
matrix & | t13, | |||
matrix & | t14, | |||
matrix & | t23, | |||
matrix & | t34, | |||
matrix & | t42, | |||
vector & | x, | |||
vector & | y, | |||
vector & | z | |||
) |
Transformation
volcoord | - volume coordinates | |
nodval | - nodal values |
20.8.2008
Definition at line 140 of file lintetrot.cpp.
References ned.
Referenced by bf_matrix(), geom_matrix(), and geom_matrix_shear().
function assembles transformation matrix from nodal to global coordinate system
nodes | - nodes of element | |
tmat | - transformation matrix |
JK, 2.10.2008
Definition at line 519 of file lintetrot.cpp.
References node::e1, node::e2, node::e3, 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_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 1496 of file lintetrot.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().
long* cncomp |
cumulative array of numbers of components of blocks
Definition at line 112 of file lintetrot.h.
Referenced by eigstrain_forces(), elem_integration(), lintetrot(), and ~lintetrot().
long intordb |
order of numerical interation on surfaces
Definition at line 104 of file lintetrot.h.
Referenced by lintetrot(), node_forces_surf(), and node_forces_surf_old().
long intordmm |
order of integration of mass matrix
Definition at line 102 of file lintetrot.h.
Referenced by lintetrot(), load_matrix(), and mass_matrix().
long** intordsm |
array of orders of integration of stiffness matrix
Definition at line 100 of file lintetrot.h.
Referenced by compute_nlstress(), dmatblock(), mechtop::give_intordsm(), inicipval(), ip_strains(), ipcoord(), lintetrot(), stiffness_matrix(), and ~lintetrot().
long napfun |
number of approximated functions on the element
Definition at line 90 of file lintetrot.h.
Referenced by mechtop::give_napfun(), lintetrot(), load_matrix(), mass_matrix(), node_forces_surf(), and node_forces_surf_old().
long nb |
number of blocks
Definition at line 108 of file lintetrot.h.
Referenced by mechtop::give_nb(), mechtop::give_nb_te(), inicipval(), lintetrot(), and ~lintetrot().
long* ncomp |
array of numbers of components of blocks
Definition at line 110 of file lintetrot.h.
Referenced by eigstrain_forces(), lintetrot(), stiffness_matrix(), and ~lintetrot().
long ndofe |
number of DOFs on the element
Definition at line 82 of file lintetrot.h.
Referenced by eigstrain_forces(), elem_integration(), mechtop::give_ndofe(), ip_strains(), lintetrot(), load_matrix(), mass_matrix(), node_forces_surf(), node_forces_surf_old(), res_eigstrain_forces(), res_incr_internal_forces(), res_internal_forces(), res_ip_strains(), res_nonloc_internal_forces(), res_stiffness_matrix(), and stiffness_matrix().
long ned |
number of edges on one element
Definition at line 92 of file lintetrot.h.
Referenced by mechtop::give_ned(), lintetrot(), and tran_side().
long** nip |
array of numbers of integration points in sets
Definition at line 106 of file lintetrot.h.
Referenced by mechtop::give_nip(), lintetrot(), and ~lintetrot().
long nne |
number of nodes on one element
Definition at line 84 of file lintetrot.h.
Referenced by eigstrain_forces(), geom_matrix(), geom_matrix_shear(), mechtop::give_nne(), inicipval(), ipcoord(), ipvolume(), lintetrot(), load_matrix(), mass_matrix(), nod_eqother_ip(), nod_strains_ip(), nod_stresses_ip(), node_forces_surf(), node_forces_surf_old(), res_eigstrain_forces(), res_incr_internal_forces(), res_internal_forces(), res_ip_strains(), res_nonloc_internal_forces(), res_stiffness_matrix(), stiffness_matrix(), and volumeip().
long nned |
number of nodes on one edge
Definition at line 94 of file lintetrot.h.
Referenced by mechtop::give_nned(), and lintetrot().
long nnsurf |
number of nodes on one surface
Definition at line 98 of file lintetrot.h.
Referenced by mechtop::give_nnsurf(), and lintetrot().
long nsurf |
number of surfaces
Definition at line 96 of file lintetrot.h.
Referenced by mechtop::give_nsurf(), and lintetrot().
stress/strain state
Definition at line 114 of file lintetrot.h.
Referenced by mechtop::give_ssst(), and lintetrot().
long tncomp |
total number of components of the strain and stress tensors
Definition at line 86 of file lintetrot.h.
Referenced by compute_eigstress(), eigstrain_forces(), elem_integration(), mechtop::give_ncomp(), mechtop::give_tncomp(), ip_strains(), lintetrot(), nod_strains_ip(), nod_stresses_ip(), res_ip_strains(), and stiffness_matrix().
long tnip |
total number of integration points on element
Definition at line 88 of file lintetrot.h.
Referenced by mechtop::give_tnip(), and lintetrot().