#include <plquadcontact.h>
Public Member Functions | |
double | approx (double xi, vector &nodval) |
void | compute_nlstress (long lcid, long eid, long ri, long ci) |
void | elem_integration (integratedquant iq, long lcid, long eid, long ri, long ci, vector &nv, vector &x, vector &y) |
void | geom_matrix (matrix &gm, vector &x, vector &y, double xi, double &jac) |
void | internal_forces (long lcid, long eid, long ri, long ci, vector &ifor, vector &x, vector &y) |
void | mainip_strains (long lcid, long eid, long ri, long ci, vector &r) |
plquadcontact (void) | |
void | res_internal_forces (long lcid, long eid, vector &ifor) |
void | res_mainip_strains (long lcid, long eid) |
void | res_mainip_stresses (long lcid, long eid) |
void | res_stiffness_matrix (long eid, matrix &sm) |
void | stiffness_matrix (long eid, long ri, long ci, matrix &sm) |
void | transf_matrix (ivector &nod, matrix &tmat) |
~plquadcontact (void) | |
Public Attributes | |
long ** | intordsm |
array containing cumulative numbers of components of stress and strain tensors | |
long | nb |
number of blocks | |
long * | ncomp |
array containing numbers of components of stress and strain tensors | |
long | ndofe |
number of DOFs on the element | |
long ** | nip |
order of integration of mass matrix | |
long | nne |
number of nodes on one element | |
strastrestate | ssst |
stress/strain state | |
long | tncomp |
total number of components of stress and strain tensors | |
long | tnip |
total number of integration points on element |
class plquadcontact is plane quadrilateral element for contact problems
basic data nne = 4 - number nodes on element ndofe = 8 - number of DOFs on element
JK, 11.6.2006
Definition at line 21 of file plquadcontact.h.
plquadcontact | ( | void | ) |
~plquadcontact | ( | void | ) |
double approx | ( | double | xi, | |
vector & | nodval | |||
) |
The function initializes basic data on elements.
eid | - element id |
JK, 11.6.2006
void plquadcontact::eleminit (long eid) {
long ii,jj;
Mt->elements[eid].nb=nb; Mt->elements[eid].intordsm = new long* [nb]; Mt->elements[eid].nip = new long* [nb];
for (ii=0;ii<nb;ii++){ Mt->elements[eid].intordsm[ii] = new long [nb]; Mt->elements[eid].nip[ii] = new long [nb]; for (jj=0;jj<nb;jj++){ Mt->elements[eid].intordsm[ii][jj]=intordsm[ii][jj]; Mt->elements[eid].nip[ii][jj]=nip[ii][jj]; } } } function approximates function defined by nodal values
xi,eta | - coordinates on element | |
nodval | - nodal values |
JK
Definition at line 88 of file plquadcontact.cpp.
References vector::a, bf_lin_1d(), f, nne, and scprd().
Referenced by elem_integration(), and stiffness_matrix().
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 |
Created by Tomas Koudelka, 11.2008
Definition at line 415 of file plquadcontact.cpp.
References mechmat::computenlstresses(), mechtop::elements, intordsm, element::ipp, Mm, Mp, Mt, and probdesc::strcomp.
Referenced by internal_forces(), and res_mainip_stresses().
void elem_integration | ( | integratedquant | iq, | |
long | lcid, | |||
long | eid, | |||
long | ri, | |||
long | ci, | |||
vector & | nv, | |||
vector & | x, | |||
vector & | y | |||
) |
The function integrates selected quantity on the selected 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 | - nodal coordinates |
JF, 29.10.2012
Definition at line 510 of file plquadcontact.cpp.
References vector::a, addv(), approx(), cmulv(), mechtop::elements, fillv(), gauss_points(), geom_matrix(), mechtop::give_elemnodes(), mechcrsec::give_thickness(), mechmat::givequantity(), intordsm, element::ipp, Mc, Mm, Mt, mtxv(), ndofe, nne, nodes, and tncomp.
Referenced by internal_forces().
The function assembles strain-displacement (geometric) matrix.
gm | - resulting geometric matrix (output) | |
x,y | - nodal coordinates | |
xi | - natural coordinate of required integration point | |
jac | - Jacobian of transformation (output) |
JF, 15.10.2012
Definition at line 322 of file plquadcontact.cpp.
References cmulm(), fillm(), and mxm().
Referenced by elem_integration(), mainip_strains(), and stiffness_matrix().
void internal_forces | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci, | |||
vector & | ifor, | |||
vector & | x, | |||
vector & | y | |||
) |
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 (output) | |
x,y | - vectors of nodal coordinates |
Created by Tomas Koudelka, 11.2008
Definition at line 445 of file plquadcontact.cpp.
References compute_nlstress(), elem_integration(), and locstress.
Referenced by res_internal_forces().
void mainip_strains | ( | long | lcid, | |
long | eid, | |||
long | ri, | |||
long | ci, | |||
vector & | r | |||
) |
The function computes strains at integration points of element. It is used in geometrically linear problems.
lcid | - load case id | |
eid | - element id | |
ri,ci | - row and column indices | |
ii | - number of block | |
r | - vector of nodal displacements |
JK, 11.6.2006
Definition at line 274 of file plquadcontact.cpp.
References mechtop::elements, fillm(), gauss_points(), geom_matrix(), mechtop::give_elemnodes(), mechtop::give_node_coord2d(), intordsm, element::ipp, Mm, Mt, mxv(), ncomp, ndofe, nne, nodes, mechmat::storestrain(), and tncomp.
Referenced by res_mainip_strains().
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 |
Created by Tomas Koudelka, 11.2008
Definition at line 470 of file plquadcontact.cpp.
References copyv(), mechtop::give_elemnodes(), mechtop::give_node_coord2d(), glvectortransf(), internal_forces(), mechtop::locsystems(), Mt, ndofe, nne, nodes, and transf_matrix().
Referenced by elem_internal_forces(), and elem_nonloc_internal_forces().
void res_mainip_strains | ( | long | lcid, | |
long | eid | |||
) |
The function computes strains at integration points. If there are defined loacal nodal coordinate systems then the corresponding components of the dislacement vector are transformed.
lcid | - load case id | |
eid | - element id |
JK, 11.6.2006
Definition at line 230 of file plquadcontact.cpp.
References vector::a, allocm(), allocv(), copyv(), destrm(), destrv(), eldispl(), lgvectortransf(), mechtop::locsystems(), mainip_strains(), Mt, ndofe, nne, nodes, and transf_matrix().
Referenced by compute_ipstrains().
void res_mainip_stresses | ( | long | lcid, | |
long | eid | |||
) |
The function computes stresses at element integration points
lcid | - load case id | |
eid | - element id |
TKo, 11.2012
Definition at line 392 of file plquadcontact.cpp.
References compute_nlstress().
Referenced by compute_ipstresses().
void res_stiffness_matrix | ( | long | eid, | |
matrix & | sm | |||
) |
Function computes stiffness matrix of one element. If it is required, nodal values are transformed to the local nodal coordinate systems.
eid | - number of element | |
sm | - stiffness matrix (output) |
Created by Tomas Koudelka, 11.2008
Definition at line 197 of file plquadcontact.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 | |||
) |
The function computes stiffness matrix of plane rectangular finite element for contact problems.
eid | - number of element | |
ri,ci | - row and column indices | |
sm | - resulting stiffness matrix |
JF, 29.10.2012
Definition at line 148 of file plquadcontact.cpp.
References vector::a, approx(), bdbjac(), mechtop::elements, fillm(), gauss_points(), geom_matrix(), mechtop::give_elemnodes(), mechtop::give_node_coord2d(), mechcrsec::give_thickness(), intordsm, element::ipp, mechmat::matstiff(), Mc, Mm, Mt, ncomp, ndofe, nne, nodes, and tncomp.
Referenced by res_stiffness_matrix().
The function assembles transformation matrix x_g = T x_l for transformation from the local nodal coordinate systems.
nodes | - array containing node numbers | |
tmat | - transformation matrix |
JK, 11.6.2006
Definition at line 114 of file plquadcontact.cpp.
References node::e1, node::e2, fillm(), matrix::m, Mt, ivector::n, mechtop::nodes, and node::transf.
Referenced by res_internal_forces(), res_mainip_strains(), and res_stiffness_matrix().
long** intordsm |
array containing cumulative numbers of components of stress and strain tensors
number of approximated functions on the element number of edges number of nodes on one edge array containing orders of numerical integration of stiffness matrix
Definition at line 69 of file plquadcontact.h.
Referenced by compute_nlstress(), elem_integration(), mechtop::give_intordsm(), mainip_strains(), plquadcontact(), stiffness_matrix(), and ~plquadcontact().
long nb |
number of blocks
Definition at line 77 of file plquadcontact.h.
Referenced by mechtop::give_nb(), mechtop::give_nb_te(), plquadcontact(), and ~plquadcontact().
long* ncomp |
array containing numbers of components of stress and strain tensors
Definition at line 59 of file plquadcontact.h.
Referenced by mainip_strains(), plquadcontact(), stiffness_matrix(), and ~plquadcontact().
long ndofe |
number of DOFs on the element
Definition at line 51 of file plquadcontact.h.
Referenced by elem_integration(), mechtop::give_ndofe(), mainip_strains(), plquadcontact(), res_internal_forces(), res_mainip_strains(), res_stiffness_matrix(), and stiffness_matrix().
long** nip |
order of integration of mass matrix
order of integration on edges array of numbers of integration points in sets
Definition at line 75 of file plquadcontact.h.
Referenced by mechtop::give_nip(), plquadcontact(), and ~plquadcontact().
long nne |
number of nodes on one element
Definition at line 53 of file plquadcontact.h.
Referenced by approx(), elem_integration(), mechtop::give_nne(), mainip_strains(), plquadcontact(), res_internal_forces(), res_mainip_strains(), res_stiffness_matrix(), and stiffness_matrix().
stress/strain state
Definition at line 79 of file plquadcontact.h.
Referenced by mechtop::give_ssst(), and plquadcontact().
long tncomp |
total number of components of stress and strain tensors
Definition at line 55 of file plquadcontact.h.
Referenced by elem_integration(), mechtop::give_ncomp(), mechtop::give_tncomp(), mainip_strains(), plquadcontact(), and stiffness_matrix().
long tnip |
total number of integration points on element
Definition at line 57 of file plquadcontact.h.
Referenced by mechtop::give_tnip(), and plquadcontact().