#include <quadhex.h>
Public Member Functions | |
double | approx (double xi, double eta, double zeta, vector &nodval) |
void | aver_strains (long lcid, long eid, long ri, long ci, vector &averstra, double &volume) |
void | bf_matrix (matrix &n, double xi, double eta, double zeta) |
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 (long ri, long ci, matrix &d, matrix &dd) |
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, double xi, double eta, double zeta, double &jac) |
void | geom_matrix_block (matrix &gm, long ri, vector &x, vector &y, vector &z, double xi, double eta, double zeta, double &jac) |
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 | intpointval2 (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 &x, vector &y, vector &z, 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) |
void | local_values (long lcid, long eid, long ri, long ci) |
void | mass_matrix (long eid, matrix &mm) |
void | nod_eqother_ip (long eid, long ri, long ci) |
void | nod_strains_comp (long lcid, long eid, double **stra) |
void | nod_strains_ip (long lcid, long eid, long ri, long ci) |
void | nod_stresses_comp (long lcid, long eid, long ri, long ci, double **stra, double **stre) |
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 | nodecoord (vector &xi, vector &eta, vector &zeta) |
void | nonloc_internal_forces (long lcid, long eid, long ri, long ci, vector &ifor, vector &x, vector &y, vector &z) |
quadhex (void) | |
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 (matrix &tran, vector &gx, vector &gy, vector &gz, long is) |
void | transf_matrix (ivector &nodes, matrix &tmat) |
~quadhex (void) | |
Public Attributes | |
long * | cncomp |
array containing cumulative numbers of components of stress and strain tensors | |
long | intordb |
array of numbers of integration points in sets | |
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 containing numbers of components of stress and strain tensors | |
long | ndofe |
number of DOFs on the element | |
long | ned |
number of edges | |
long ** | nip |
array of numbers of integration points on surface | |
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 |
number of components of the strain and stress tensors | |
long | tnip |
total number of integration points on element |
class quadhex defines hexahedral finite element with tri-quadratic approximation functions
JK
Definition at line 16 of file quadhex.h.
quadhex | ( | void | ) |
~quadhex | ( | void | ) |
double approx | ( | double | xi, | |
double | eta, | |||
double | zeta, | |||
vector & | nodval | |||
) |
function approximates function defined by nodal values
xi,eta,zeta | - natural coordinates | |
nodval | - nodal values |
JK, 20.8.2001
Definition at line 98 of file quadhex.cpp.
References vector::a, bf_quad_hex_3d(), f, nne, and scprd().
Referenced by inicipval(), interpolate(), intpointval(), ipcoord(), and mass_matrix().
void aver_strains | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci, | |||
vector & | averstra, | |||
double & | volume | |||
) |
Definition at line 1988 of file quadhex.cpp.
References vector::a, allocv(), destrv(), mechtop::elements, gauss_points(), mechtop::give_node_coord3d(), mechmat::givestress(), intordsm, element::ipp, jac_3d(), Mm, Mt, vector::n, nb, ncomp, and nne.
Referenced by mechmat::compute_averstrains().
void bf_matrix | ( | matrix & | n, | |
double | xi, | |||
double | eta, | |||
double | zeta | |||
) |
function assembles matrix of base functions
n | - matrix of base functions | |
xi,eta,zeta | - coordinates |
JK, 16.8.2001
Definition at line 118 of file quadhex.cpp.
References vector::a, bf_quad_hex_3d(), fillm(), and nne.
Referenced by load_matrix(), mass_matrix(), and node_forces_surf().
void compute_eigstress | ( | 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 |
JK, 27.11.2006
Definition at line 1274 of file quadhex.cpp.
References mechtop::elements, mechmat::giveeigstrain(), intordsm, 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 |
JK, 27.11.2006
Definition at line 1157 of file quadhex.cpp.
References mechmat::computenlstresses(), mechtop::elements, intordsm, 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
lcid | - number of load case | |
eid | - element id | |
ri,ci | - row and column indices |
TKo 7.2008
Definition at line 1186 of file quadhex.cpp.
References mechmat::computenlstressesincr(), mechtop::elements, intordsm, 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 |
JK, 27.11.2006
Definition at line 1247 of file quadhex.cpp.
References mechmat::compnonloc_nlstresses(), mechtop::elements, intordsm, 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, | |||
vector & | x, | |||
vector & | y, | |||
vector & | z | |||
) |
function computes contributions from eigenstrains to the right hand side
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 |
JK, 17.8.2004 TKo 7.2008
Definition at line 1006 of file quadhex.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 |
JK, 27.11.2006 TKo 7.2008
Definition at line 1312 of file quadhex.cpp.
References vector::a, addv(), allocv(), cmulv(), cncomp, destrv(), mechtop::elements, fillv(), gauss_points(), geom_matrix(), mechmat::givequantity(), intordsm, element::ipp, Mm, Mt, mtxv(), ndofe, and tncomp.
Referenced by eigstrain_forces(), incr_internal_forces(), internal_forces(), and nonloc_internal_forces().
void geom_matrix | ( | matrix & | gm, | |
vector & | x, | |||
vector & | y, | |||
vector & | z, | |||
double | xi, | |||
double | eta, | |||
double | zeta, | |||
double & | jac | |||
) |
function computes strain-displacement (geometric) matrix
gm | - geometric matrix | |
x,y,z | - vectors containing element node coordinates | |
xi,eta,zeta | - coordinates | |
jac | - Jacobian |
JK, 16.8.2001
Definition at line 145 of file quadhex.cpp.
References vector::a, derivatives_3d(), dx_bf_quad_hex_3d(), dy_bf_quad_hex_3d(), dz_bf_quad_hex_3d(), fillm(), and nne.
Referenced by elem_integration(), ip_strains(), hexahedc::lower_cap_coup_matrix(), hexahedc::lower_cond_coup_matrix(), nod_strains_comp(), stiffness_matrix(), hexahedc::upper_cap_coup_matrix(), hexahedc::upper_cond_coup_matrix(), and hexahedc::upper_cond_coup_vector().
void geom_matrix_block | ( | matrix & | gm, | |
long | ri, | |||
vector & | x, | |||
vector & | y, | |||
vector & | z, | |||
double | xi, | |||
double | eta, | |||
double | zeta, | |||
double & | jac | |||
) |
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 981 of file quadhex.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 1412 of file quadhex.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, | |||
vector & | x, | |||
vector & | y, | |||
vector & | z | |||
) |
function computes internal forces in the case of geometrical linear computation
lcid | - number of load case | |
eid | - element id | |
ri,ci | - row and column indices | |
ifor | - vector of internal forces | |
x,y,z | - vectors of nodal coordinates |
JK, 28.7.2001 TKo 7.2008
Definition at line 929 of file quadhex.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 quadratic approximation functions are used
eid | - element id | |
nodval | - nodal values | |
ipval | - value at integration points |
17.8.2004, JK
Definition at line 1901 of file quadhex.cpp.
References vector::a, allocv(), approx(), destrv(), gauss_points(), intordsm, and nb.
Referenced by elem_intpointval().
function interpolates the nodal values to the integration points on the element linear approximation functions are used
eid | - element id | |
nodval | - nodal values | |
ipval | - value at integration points |
17.8.2004, JK
Definition at line 1946 of file quadhex.cpp.
References vector::a, allocv(), linhex::approx(), destrv(), gauss_points(), intordsm, Lhex, nb, and linhex::nne.
Referenced by elem_intpointval2().
void ip_elast_stresses | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci | |||
) |
function computes stresses at integration points of element stresses are computed from strains with the help of elastic stiffness
lcid | - load case id | |
eid | - element id | |
ri | - row index | |
ci | - column index |
JK, 10.5.2002
Definition at line 721 of file quadhex.cpp.
References vector::a, allocv(), destrv(), mechtop::elements, gauss_points(), mechmat::givestrain(), intordsm, element::ipp, mechmat::matstiff(), Mm, Mt, mxv(), ndofe, mechmat::storestress(), and tncomp.
void ip_strains | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci, | |||
vector & | x, | |||
vector & | y, | |||
vector & | z, | |||
vector & | r | |||
) |
function computes strains at integration points of element
lcid | - load case id | |
eid | - element id | |
ri | - row index | |
ci | - column index | |
x,y,z | - vectors of nodal coordinates | |
r | - vector of nodal displacements |
JK, 10.5.2002, modified 27.11.2006
Definition at line 484 of file quadhex.cpp.
References vector::a, allocm(), allocv(), cncomp, destrm(), destrv(), mechtop::elements, gauss_points(), geom_matrix(), intordsm, element::ipp, Mm, Mt, mxv(), nb, ncomp, ndofe, nne, nodes, and mechmat::storestrain().
Referenced by res_ip_strains().
void ip_stresses | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci | |||
) |
function computes stresses at integration points of element stresses are computed by material models
lcid | - load case id | |
eid | - element id | |
ri | - row index | |
ci | - column index |
10.5.2002, JK
Definition at line 690 of file quadhex.cpp.
References mechmat::computenlstresses(), mechtop::elements, intordsm, 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
eid | - element id | |
ipp | - integration point pointer | |
ri | - row index | |
ci | - column index | |
coord | - vector of coordinates |
10.1.2002
Definition at line 1361 of file quadhex.cpp.
References vector::a, approx(), mechtop::elements, gauss_points(), 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 1532 of file quadhex.cpp.
References vector::a, allocv(), destrv(), mechtop::elements, gauss_points(), mechtop::give_node_coord3d(), intordsm, element::ipp, jac_3d(), Mm, Mt, nb, 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 |
JK, 16.8.2001
Definition at line 375 of file quadhex.cpp.
References matrix::a, vector::a, bf_matrix(), fillm(), gauss_points(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), intordmm, jac_3d(), matrix::m, Mt, matrix::n, napfun, ndofe, nne, nnj(), and nodes.
Referenced by loadmat(), and 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.
lcid | - number of load case | |
eid | - element id | |
ri,ci | - row and column indices |
TKo 7.2008
Definition at line 1218 of file quadhex.cpp.
References mechmat::computenlstresses(), mechtop::elements, intordsm, element::ipp, Mm, Mp, Mt, and probdesc::strcomp.
Referenced by elem_local_values().
void mass_matrix | ( | long | eid, | |
matrix & | mm | |||
) |
function computes mass matrix
eid | - number of element | |
mm | - mass matrix |
JK, 16.8.2001
Definition at line 303 of file quadhex.cpp.
References matrix::a, approx(), bf_matrix(), fillm(), gauss_points(), mechcrsec::give_density(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), intordmm, jac_3d(), matrix::m, Mc, Mt, matrix::n, napfun, ndofe, nne, nnj(), and nodes.
Referenced by massmat(), and res_mass_matrix().
void nod_eqother_ip | ( | long | eid, | |
long | ri, | |||
long | ci | |||
) |
function computes other values in nodes of element
eid | - element id | |
ri,ci | - row and column indices |
JK, 24.10.2005
Definition at line 884 of file quadhex.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_quadhex(), and node::storeother().
Referenced by compute_nodeothers().
void nod_strains_comp | ( | long | lcid, | |
long | eid, | |||
double ** | stra | |||
) |
function computes strains in nodes of element
lcid | - load case id | |
eid | - element id | |
stra | - array for strain components |
stra[i][j] - the j-th strain component at the i-th node
JK, 10.5.2002
Definition at line 570 of file quadhex.cpp.
References allocm(), allocv(), copyv(), destrm(), destrv(), eldispl(), geom_matrix(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), lgvectortransf(), mechtop::locsystems(), Mt, mxv(), ndofe, nne, nodcoord_quadhex(), nodes, tncomp, and transf_matrix().
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 | |
ri,ci | - row and column indices |
JK, 27.9.2005
Definition at line 534 of file quadhex.cpp.
References mechtop::elements, mechtop::give_elemnodes(), mechmat::givestrain(), intordsm, element::ipp, Mm, Mt, nne, nod, mechtop::nodes, nodip_quadhex(), node::storestrain(), and tncomp.
Referenced by compute_nodestrains().
void nod_stresses_comp | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci, | |||
double ** | stra, | |||
double ** | stre | |||
) |
function computes stresses in nodes
lcid | - load case id | |
eid | - element id | |
ri,ci | - row and column indices | |
stra | - | |
stre | - |
10.5.2002
Definition at line 802 of file quadhex.cpp.
References mechtop::elements, element::ipp, mechmat::matstiff(), Mm, Mt, mxv(), nne, and tncomp.
void nod_stresses_ip | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci | |||
) |
function computes stresses at nodes of element
lcid | - load case id | |
eid | - element id | |
ri,ci | - row and column indices |
JK, 27.9.2005
Definition at line 766 of file quadhex.cpp.
References mechtop::elements, mechtop::give_elemnodes(), mechmat::givestress(), intordsm, element::ipp, Mm, Mt, nne, nod, mechtop::nodes, nodip_quadhex(), node::storestress(), and tncomp.
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
eid | - element id | |
ri,ci | - row and column indices | |
nfor | - vector of presure | |
eis | - surface id |
4.2002, PF
Definition at line 1584 of file quadhex.cpp.
References addv(), bf_matrix(), fillm(), fillv(), gauss_points(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), glvectortransfblock(), intordb, jac2d3d(), lgvectortransfblock(), Mt, mxv(), napfun, ndofe, nne, nnj(), nnsurf, nodes, and tran_mat().
Referenced by loadel::surfaceload().
void nonloc_internal_forces | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci, | |||
vector & | ifor, | |||
vector & | x, | |||
vector & | y, | |||
vector & | z | |||
) |
function computes nonlocal internal forces
lcid | - number of load case | |
eid | - element id | |
ri,ci | - row and column indices | |
ifor | - vector of internal forces | |
x,y,z | - nodal coordinates |
JK, 28.7.2001 TKo 7.2008
Definition at line 955 of file quadhex.cpp.
References compute_nonloc_nlstress(), elem_integration(), and nonlocstress.
Referenced by res_nonloc_internal_forces().
void res_eigstrain_forces | ( | long | lcid, | |
long | eid, | |||
vector & | nfor | |||
) |
function computes contributions from eigenstrains to the right hand side
eid | - element id | |
ifor | - vector of internal forces |
JK, 17.8.2004
Definition at line 1123 of file quadhex.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 | |||
) |
Definition at line 1090 of file quadhex.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 (from correct stresses)
lcid | - number of load case | |
eid | - element id | |
ifor | - vector of internal forces |
JK, 24.9.2005
Definition at line 1031 of file quadhex.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
lcid | - load case id | |
eid | - element id |
JK, modified 23.11.2006
Definition at line 446 of file quadhex.cpp.
References allocm(), allocv(), copyv(), destrm(), destrv(), eldispl(), mechtop::give_elemnodes(), mechtop::give_node_coord3d(), 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 stresses at integration points of element
lcid | - load case id | |
eid | - element id |
JK, 10.5.2002
Definition at line 674 of file quadhex.cpp.
References ip_stresses().
Referenced by compute_ipstresses().
void res_load_matrix | ( | long | eid, | |
matrix & | lm | |||
) |
function computes load matrix
eid | - number of element | |
lm | - load matrix |
JK, 16.8.2001
Definition at line 417 of file quadhex.cpp.
References mechtop::give_elemnodes(), glmatrixtransf(), load_matrix(), mechtop::locsystems(), Mt, ndofe, nne, nodes, and transf_matrix().
void res_mass_matrix | ( | long | eid, | |
matrix & | mm | |||
) |
function computes mass matrix
eid | - number of element | |
mm | - mass matrix |
JK, 16.8.2001
Definition at line 349 of file quadhex.cpp.
References mechtop::give_elemnodes(), 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 internal forces for nonlocal models
lcid | - load case id | |
eid | - element id | |
ifor | - vector of internal forces |
TKo, 7.2008
Definition at line 1065 of file quadhex.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 computes stiffness matrix of one element
eid | - number of element | |
sm | - stiffness matrix |
JK, 16.8.2001
Definition at line 277 of file quadhex.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 |
JK, 16.8.2001
Definition at line 227 of file quadhex.cpp.
References vector::a, allocv(), bdbjac(), destrv(), mechtop::elements, fillm(), gauss_points(), geom_matrix(), mechtop::give_node_coord3d(), intordsm, 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 | |||
) |
Definition at line 619 of file quadhex.cpp.
References allocv(), destrv(), enodes, aepoints::give_aepcoord(), aepoints::give_naep(), aepoints::give_ncomp(), aepoints::give_sid(), intpts, Mm, Mp, nowhere, aepoints::storevalues(), mechmat::stra, probdesc::strainaver, aepoints::tape, and userdefined.
void stresses | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci | |||
) |
Definition at line 831 of file quadhex.cpp.
References allocv(), destrv(), enodes, aepoints::give_aepcoord(), aepoints::give_naep(), aepoints::give_ncomp(), aepoints::give_sid(), intpts, Mm, Mp, nowhere, aepoints::storevalues(), mechmat::stre, probdesc::stressaver, aepoints::tape, and userdefined.
Referenced by computestresses().
function computes transformation matrix on surface
x,y | - local coordinate xL=adge12 | |
tran | - tranformation metrix to GCS | |
gx,gy,gz | - vector of node global coordinate | |
is | - surface id |
4.2002, PF
Definition at line 1799 of file quadhex.cpp.
Referenced by node_forces_surf().
function assembles transformation matrix from local nodal coordinate system to the global coordinate system x_g = T x_l
nodes | - nodes of element | |
tmat | - transformation matrix |
JK
Definition at line 189 of file quadhex.cpp.
References node::e1, node::e2, node::e3, fillm(), matrix::m, Mt, ivector::n, mechtop::nodes, and node::transf.
Referenced by nod_strains_comp(), 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().
long* cncomp |
array containing cumulative numbers of components of stress and strain tensors
Definition at line 93 of file quadhex.h.
Referenced by elem_integration(), ip_strains(), quadhex(), and ~quadhex().
long intordb |
array of numbers of integration points in sets
Definition at line 109 of file quadhex.h.
Referenced by node_forces_surf(), and quadhex().
long intordmm |
order of integration of mass matrix
Definition at line 107 of file quadhex.h.
Referenced by load_matrix(), mass_matrix(), and quadhex().
long** intordsm |
array of orders of integration of stiffness matrix
Definition at line 105 of file quadhex.h.
Referenced by aver_strains(), compute_eigstress(), compute_nlstress(), compute_nlstressincr(), compute_nonloc_nlstress(), elem_integration(), export_gid_gauss_pt(), mechtop::give_intordsm(), inicipval(), intpointval(), intpointval2(), ip_elast_stresses(), ip_strains(), ip_stresses(), ipcoord(), ipvolume(), local_values(), nod_eqother_ip(), nod_strains_ip(), nod_stresses_ip(), quadhex(), stiffness_matrix(), and ~quadhex().
long napfun |
number of approximated functions on the element
Definition at line 95 of file quadhex.h.
Referenced by mechtop::give_napfun(), load_matrix(), mass_matrix(), node_forces_surf(), and quadhex().
long nb |
number of blocks
Definition at line 113 of file quadhex.h.
Referenced by aver_strains(), export_gid_gauss_pt(), mechtop::give_nb(), mechtop::give_nb_te(), hexahedc::hexahedc(), inicipval(), intpointval(), intpointval2(), ip_strains(), ipvolume(), quadhex(), and ~quadhex().
long* ncomp |
array containing numbers of components of stress and strain tensors
Definition at line 91 of file quadhex.h.
Referenced by aver_strains(), hexahedc::hexahedc(), ip_strains(), quadhex(), and ~quadhex().
long ndofe |
number of DOFs on the element
Definition at line 83 of file quadhex.h.
Referenced by elem_integration(), mechtop::give_ndofe(), hexahedc::hexahedc(), ip_elast_stresses(), ip_strains(), load_matrix(), mass_matrix(), nod_strains_comp(), node_forces_surf(), quadhex(), res_eigstrain_forces(), res_incr_internal_forces(), res_internal_forces(), res_ip_strains(), res_load_matrix(), res_mass_matrix(), res_nonloc_internal_forces(), res_stiffness_matrix(), and stiffness_matrix().
long ned |
number of edges
Definition at line 97 of file quadhex.h.
Referenced by mechtop::give_ned(), and quadhex().
long** nip |
array of numbers of integration points on surface
Definition at line 111 of file quadhex.h.
Referenced by mechtop::give_nip(), quadhex(), and ~quadhex().
long nne |
number of nodes on one element
Definition at line 85 of file quadhex.h.
Referenced by approx(), aver_strains(), bf_matrix(), geom_matrix(), mechtop::give_nne(), hexahedc::hexahedc(), inicipval(), ip_strains(), ipcoord(), ipvolume(), load_matrix(), mass_matrix(), nod_eqother_ip(), nod_strains_comp(), nod_strains_ip(), nod_stresses_comp(), nod_stresses_ip(), node_forces_surf(), quadhex(), res_eigstrain_forces(), res_incr_internal_forces(), res_internal_forces(), res_ip_strains(), res_load_matrix(), res_mass_matrix(), res_nonloc_internal_forces(), res_stiffness_matrix(), and stiffness_matrix().
long nned |
number of nodes on one edge
Definition at line 99 of file quadhex.h.
Referenced by mechtop::give_nned(), and quadhex().
long nnsurf |
number of nodes on one surface
Definition at line 103 of file quadhex.h.
Referenced by mechtop::give_nnsurf(), node_forces_surf(), and quadhex().
long nsurf |
number of surfaces
Definition at line 101 of file quadhex.h.
Referenced by mechtop::give_nsurf(), and quadhex().
stress/strain state
Definition at line 115 of file quadhex.h.
Referenced by mechtop::give_ssst(), and quadhex().
long tncomp |
number of components of the strain and stress tensors
Definition at line 87 of file quadhex.h.
Referenced by compute_eigstress(), elem_integration(), mechtop::give_ncomp(), mechtop::give_tncomp(), hexahedc::hexahedc(), ip_elast_stresses(), nod_strains_comp(), nod_strains_ip(), nod_stresses_comp(), nod_stresses_ip(), quadhex(), and stiffness_matrix().
long tnip |
total number of integration points on element
Definition at line 89 of file quadhex.h.
Referenced by mechtop::give_tnip(), and quadhex().