FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
m_contact_lib Module Reference

This module provide functions of contact stiffness calculation. More...

Data Types

type  tcontactstate
 This structure records contact status. More...
 

Functions/Subroutines

subroutine contact_state_init (cstate)
 Initializer. More...
 
subroutine contact_state_copy (cstate1, cstate2)
 Copy. More...
 
subroutine print_contact_state (fnum, cstate)
 Print out contact state. More...
 
subroutine contact2mpcval (cstate, etype, nnode, mpcval)
 Transfer contact condition int mpc bundary conditions. More...
 
subroutine contact2stiff (flag, cstate, etype, nnode, ele, mu, mut, fcoeff, symm, stiff, force)
 This subroutine calculate contact stiff matrix and contact force. More...
 
subroutine tied2stiff (flag, cstate, etype, nnode, mu, mut, stiff, force)
 This subroutine calculate contact stiff matrix and contact force. More...
 
subroutine getmetrictensor (pos, etype, ele, tensor)
 This subroutine calculate the metric tensor of a elemental surface. More...
 
subroutine dispincrematrix (pos, etype, nnode, ele, tangent, tensor, matrix)
 This subroutine calculate the relation between global disp and displacement along natural coordinate of master surface supposing penetration is small. More...
 
subroutine project_point2element (xyz, etype, nn, elemt, reflen, cstate, isin, distclr, ctpos, localclr)
 This subroutine find the projection of a slave point onto master surface. More...
 
subroutine project_point2solidelement (xyz, etype, nn, elemt, reflen, cstate, isin, distclr, ctpos, localclr)
 This subroutine find the projection of a slave point onto master surface. More...
 
subroutine update_tangentforce (etype, nn, elemt0, elemt, cstate)
 This subroutine find the projection of a slave point onto master surface. More...
 

Variables

integer, parameter contactunknown = -1
 
integer, parameter contactfree = -1
 contact state definition More...
 
integer, parameter contactstick = 1
 
integer, parameter contactslip = 2
 
integer, parameter contacttied = 1
 contact type or algorithm definition More...
 
integer, parameter contactglued = 2
 
integer, parameter contactsslid = 3
 
integer, parameter contactfslid = 4
 

Detailed Description

This module provide functions of contact stiffness calculation.

Function/Subroutine Documentation

◆ contact2mpcval()

subroutine m_contact_lib::contact2mpcval ( type(tcontactstate), intent(in)  cstate,
integer, intent(in)  etype,
integer, intent(in)  nnode,
real(kind=kreal), dimension(nnode*3 + 4), intent(out)  mpcval 
)

Transfer contact condition int mpc bundary conditions.

Parameters
[in]cstatecontact state
[in]etypetype of contacting surface
[in]nnodenumber of elemental nodes
[out]mpcvalMPC constraint

Definition at line 73 of file contact_lib.f90.

Here is the call graph for this function:

◆ contact2stiff()

subroutine m_contact_lib::contact2stiff ( integer, intent(in)  flag,
type(tcontactstate), intent(in)  cstate,
integer, intent(in)  etype,
integer, intent(in)  nnode,
real(kind=kreal), dimension(3,nnode), intent(in)  ele,
real(kind=kreal), intent(in)  mu,
real(kind=kreal), intent(in)  mut,
real(kind=kreal), intent(in)  fcoeff,
logical, intent(in)  symm,
real(kind=kreal), dimension(:,:), intent(out)  stiff,
real(kind=kreal), dimension(:), intent(out)  force 
)

This subroutine calculate contact stiff matrix and contact force.

Parameters
[in]flagsmall slid or finite slide
[in]cstatecontact state
[in]etypetype of contacting surface
[in]nnodenumber of elemental nodes
[in]elecoord of surface element
[in]mupenalty
[in]mutpenalty along tangent
[in]fcoefffriction coefficient
[in]symmsymmtricalize
[out]stiffcontact stiffness
[out]forcecontact force

Definition at line 94 of file contact_lib.f90.

Here is the call graph for this function:

◆ contact_state_copy()

subroutine m_contact_lib::contact_state_copy ( type(tcontactstate), intent(in)  cstate1,
type(tcontactstate), intent(inout)  cstate2 
)

Copy.

Parameters
[in]cstate1contact state
[in,out]cstate2contact state

Definition at line 55 of file contact_lib.f90.

◆ contact_state_init()

subroutine m_contact_lib::contact_state_init ( type(tcontactstate), intent(inout)  cstate)

Initializer.

Parameters
[in,out]cstatecontact state

Definition at line 48 of file contact_lib.f90.

Here is the caller graph for this function:

◆ dispincrematrix()

subroutine m_contact_lib::dispincrematrix ( real(kind=kreal), dimension(2), intent(in)  pos,
integer, intent(in)  etype,
integer, intent(in)  nnode,
real(kind=kreal), dimension(3,nnode), intent(in)  ele,
real(kind=kreal), dimension(3,2), intent(out)  tangent,
real(kind=kreal), dimension(2,2), intent(out)  tensor,
real(kind=kreal), dimension(2,nnode*3+3), intent(out)  matrix 
)

This subroutine calculate the relation between global disp and displacement along natural coordinate of master surface supposing penetration is small.

Parameters
[in]poscurrent position(local coordinate)
[in]etypesurface element type
[in]nnodenumber of nodes in surface
[in]eleelemental coordinates
[out]tangenttangent basis
[out]tensormetric tensor
[out]matrixrelation between local and global disp increment

Definition at line 216 of file contact_lib.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getmetrictensor()

subroutine m_contact_lib::getmetrictensor ( real(kind=kreal), dimension(2), intent(in)  pos,
integer, intent(in)  etype,
real(kind=kreal), dimension(:,:), intent(in)  ele,
real(kind=kreal), dimension(2,2), intent(out)  tensor 
)

This subroutine calculate the metric tensor of a elemental surface.

Parameters
[in]poscurrent position(local coordinate)
[in]etypesurface element type
[in]eleelemental coordinates
[out]tensormetric tensor

Definition at line 198 of file contact_lib.f90.

Here is the call graph for this function:

◆ print_contact_state()

subroutine m_contact_lib::print_contact_state ( integer, intent(in)  fnum,
type(tcontactstate), intent(in)  cstate 
)

Print out contact state.

Parameters
[in]fnumfile number
[in]cstatecontact state

Definition at line 62 of file contact_lib.f90.

◆ project_point2element()

subroutine m_contact_lib::project_point2element ( real(kind=kreal), dimension(3), intent(in)  xyz,
integer, intent(in)  etype,
integer, intent(in)  nn,
real(kind=kreal), dimension(:,:), intent(in)  elemt,
real(kind=kreal), intent(in)  reflen,
type(tcontactstate), intent(inout)  cstate,
logical, intent(out)  isin,
real(kind=kreal), intent(in)  distclr,
real(kind=kreal), dimension(2), optional  ctpos,
real(kind=kreal), optional  localclr 
)

This subroutine find the projection of a slave point onto master surface.

Parameters
[in]xyzCoordinates of a spacial point, whose projecting point is to be computed
[in]etypesurface element type
[in]nnnumber of elemental nodes
[in]elemtnodes coordinates of surface element
[in]reflenreference length of surface element
[in,out]cstateRecorde of contact information
[out]isinin contact or not
[in]distclrclearance of contact distance
ctposcurr contact position( natural coord )
localclrclearance of contact local coord

Definition at line 257 of file contact_lib.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ project_point2solidelement()

subroutine m_contact_lib::project_point2solidelement ( real(kind=kreal), dimension(3), intent(in)  xyz,
integer, intent(in)  etype,
integer, intent(in)  nn,
real(kind=kreal), dimension(:,:), intent(in)  elemt,
real(kind=kreal), intent(in)  reflen,
type(tcontactstate), intent(inout)  cstate,
logical, intent(out)  isin,
real(kind=kreal), intent(in)  distclr,
real(kind=kreal), dimension(3), optional  ctpos,
real(kind=kreal), optional  localclr 
)

This subroutine find the projection of a slave point onto master surface.

Parameters
[in]xyzCoordinates of a spacial point, whose projecting point is to be computed
[in]etypesurface element type
[in]nnnumber of elemental nodes
[in]elemtnodes coordinates of surface element
[in]reflenreference length of surface element
[in,out]cstateRecorde of contact information
[out]isinin contact or not
[in]distclrclearance of contact distance
ctposcurr contact position( natural coord )
localclrclearance of contact local coord

Definition at line 365 of file contact_lib.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ tied2stiff()

subroutine m_contact_lib::tied2stiff ( integer, intent(in)  flag,
type(tcontactstate), intent(in)  cstate,
integer, intent(in)  etype,
integer, intent(in)  nnode,
real(kind=kreal), intent(in)  mu,
real(kind=kreal), intent(in)  mut,
real(kind=kreal), dimension(:,:), intent(out)  stiff,
real(kind=kreal), dimension(:), intent(out)  force 
)

This subroutine calculate contact stiff matrix and contact force.

Parameters
[in]flagsmall slid or finite slide
[in]cstatecontact state
[in]etypetype of contacting surface
[in]nnodenumber of elemental nodes
[in]mupenalty
[in]mutpenalty along tangent
[out]stiffcontact stiffness
[out]forcecontact force

Definition at line 166 of file contact_lib.f90.

Here is the call graph for this function:

◆ update_tangentforce()

subroutine m_contact_lib::update_tangentforce ( integer, intent(in)  etype,
integer, intent(in)  nn,
real(kind=kreal), dimension(3,nn), intent(in)  elemt0,
real(kind=kreal), dimension(3,nn), intent(in)  elemt,
type(tcontactstate), intent(inout)  cstate 
)

This subroutine find the projection of a slave point onto master surface.

Parameters
[in]etypesurface element type
[in]nnnumber of elemental nodes
[in]elemt0nodes coordinates of surface element at t
[in]elemtnodes coordinates of surface element at t+dt
[in,out]cstateRecorde of contact information

Definition at line 452 of file contact_lib.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ contactfree

integer, parameter m_contact_lib::contactfree = -1

contact state definition

Definition at line 16 of file contact_lib.f90.

◆ contactfslid

integer, parameter m_contact_lib::contactfslid = 4

Definition at line 24 of file contact_lib.f90.

◆ contactglued

integer, parameter m_contact_lib::contactglued = 2

Definition at line 22 of file contact_lib.f90.

◆ contactslip

integer, parameter m_contact_lib::contactslip = 2

Definition at line 18 of file contact_lib.f90.

◆ contactsslid

integer, parameter m_contact_lib::contactsslid = 3

Definition at line 23 of file contact_lib.f90.

◆ contactstick

integer, parameter m_contact_lib::contactstick = 1

Definition at line 17 of file contact_lib.f90.

◆ contacttied

integer, parameter m_contact_lib::contacttied = 1

contact type or algorithm definition

Definition at line 21 of file contact_lib.f90.

◆ contactunknown

integer, parameter m_contact_lib::contactunknown = -1

Definition at line 14 of file contact_lib.f90.