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

This module provides functions to calculate contact stiff matrix. More...

Functions/Subroutines

subroutine print_contatct_pair (file, pair)
 Write out the contact definition read from mesh file. More...
 
subroutine fstr_set_contact_penalty (maxv)
 
logical function fstr_is_contact_active ()
 
subroutine fstr_set_contact_active (a)
 
logical function fstr_is_contact_conv (ctAlgo, infoCTChange, hecMESH)
 
logical function fstr_is_matrixstructure_changed (infoCTChange)
 
subroutine fstr_contact2mpc (contacts, mpcs)
 Contact states to equation conditions. More...
 
subroutine fstr_del_contactmpc (mpcs)
 Delete mpcs derived from contact conditions. More...
 
subroutine fstr_write_mpc (file, mpcs)
 Print out mpc conditions. More...
 
subroutine fstr_scan_contact_state (cstep, sub_step, cont_step, dt, ctAlgo, hecMESH, fstrSOLID, infoCTChange, B)
 Scanning contact state. More...
 
subroutine remove_duplication_tiedcontact (cstep, hecMESH, fstrSOLID, infoCTChange)
 Scanning contact state. More...
 
subroutine fstr_scan_contact_state_exp (cstep, hecMESH, fstrSOLID, infoCTChange)
 Scanning contact state. More...
 
subroutine fstr_update_contact0 (hecMESH, fstrSOLID, B)
 Update lagrangian multiplier. More...
 
subroutine fstr_update_contact_multiplier (hecMESH, fstrSOLID, ctchanged)
 Update lagrangian multiplier. More...
 
subroutine fstr_ass_load_contactalag (hecMESH, fstrSOLID, B)
 Update lagrangian multiplier. More...
 
subroutine fstr_update_contact_tangentforce (fstrSOLID)
 Update tangent force. More...
 
subroutine fstr_contactbc (cstep, iter, hecMESH, hecMAT, fstrSOLID)
 Introduce contact stiff into global stiff matrix or mpc conditions into hecMESH. More...
 
subroutine initialize_contact_output_vectors (fstrSOLID, hecMAT)
 
subroutine initialize_embed_vectors (fstrSOLID, hecMAT)
 
subroutine setup_contact_elesurf_for_area (cstep, hecMESH, fstrSOLID)
 
subroutine calc_contact_area (hecMESH, fstrSOLID, flag)
 
subroutine calc_nodalarea_surfelement (etype, nn, ecoord, sid, vect)
 
subroutine fstr_setup_parancon_contactvalue (hecMESH, ndof, vec, vtype)
 

Variables

integer(kind=kint), save n_contact_mpc
 
real(kind=kreal), save mu =1.d10
 penalty, default value More...
 
real(kind=kreal), save mut =1.d6
 penalty along tangent direction More...
 
real(kind=kreal), save cdotp =1.d3
 mu=cdotp*maxval More...
 
real(kind=kreal), save cgn =1.d-5
 convergent condition of penetration More...
 
real(kind=kreal), save cgt =1.d-3
 convergent condition of relative tangent disp More...
 
real(kind=kreal), dimension(2), save gnt
 1:current average penetration; 2:current relative tangent displacement More...
 
real(kind=kreal), dimension(2), save bakgnt
 1:current average penetration; 2:current relative tangent displacement! More...
 

Detailed Description

This module provides functions to calculate contact stiff matrix.

Function/Subroutine Documentation

◆ calc_contact_area()

subroutine mcontact::calc_contact_area ( type( hecmwst_local_mesh ), intent(in)  hecMESH,
type(fstr_solid), intent(inout)  fstrSOLID,
integer(kind=kint), intent(in)  flag 
)
Parameters
[in]hecmeshtype mesh
[in,out]fstrsolidtype fstr_solid
[in]flagset 1 if called in NR iteration

Definition at line 675 of file fstr_contact.f90.

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

◆ calc_nodalarea_surfelement()

subroutine mcontact::calc_nodalarea_surfelement ( integer(kind=kint), intent(in)  etype,
integer(kind=kint), intent(in)  nn,
real(kind=kreal), dimension(:,:), intent(in)  ecoord,
integer(kind=kint), intent(in)  sid,
real(kind=kreal), dimension(:), intent(out)  vect 
)

Definition at line 726 of file fstr_contact.f90.

Here is the caller graph for this function:

◆ fstr_ass_load_contactalag()

subroutine mcontact::fstr_ass_load_contactalag ( type( hecmwst_local_mesh ), intent(in)  hecMESH,
type(fstr_solid), intent(inout)  fstrSOLID,
real(kind=kreal), dimension(:), intent(inout)  B 
)

Update lagrangian multiplier.

Parameters
[in,out]bnodal force residual

Definition at line 447 of file fstr_contact.f90.

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

◆ fstr_contact2mpc()

subroutine mcontact::fstr_contact2mpc ( type( tcontact ), dimension(:), intent(in)  contacts,
type( hecmwst_mpc ), intent(inout)  mpcs 
)

Contact states to equation conditions.

Parameters
[in]contactscurrent contact state
[in,out]mpcsto who mpc be appended

Definition at line 165 of file fstr_contact.f90.

◆ fstr_contactbc()

subroutine mcontact::fstr_contactbc ( integer(kind=kint)  cstep,
integer(kind=kint), intent(in)  iter,
type (hecmwst_local_mesh), intent(inout)  hecMESH,
type (hecmwst_matrix), intent(inout)  hecMAT,
type(fstr_solid), intent(inout)  fstrSOLID 
)

Introduce contact stiff into global stiff matrix or mpc conditions into hecMESH.

Parameters
cstepcurrent loading step
[in]iterNR iterations
[in,out]hecmeshtype mesh
[in,out]hecmattype matrix
[in,out]fstrsolidtype fstr_solid

Definition at line 480 of file fstr_contact.f90.

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

◆ fstr_del_contactmpc()

subroutine mcontact::fstr_del_contactmpc ( type( hecmwst_mpc ), intent(inout)  mpcs)

Delete mpcs derived from contact conditions.

Parameters
[in,out]mpcsmpcs to be modified

Definition at line 187 of file fstr_contact.f90.

Here is the call graph for this function:

◆ fstr_is_contact_active()

logical function mcontact::fstr_is_contact_active

Definition at line 52 of file fstr_contact.f90.

Here is the caller graph for this function:

◆ fstr_is_contact_conv()

logical function mcontact::fstr_is_contact_conv ( integer(kind=kint), intent(in)  ctAlgo,
type (fstr_info_contactchange), intent(in)  infoCTChange,
type (hecmwst_local_mesh), intent(in)  hecMESH 
)
Parameters
[in]ctalgocontact analysis algorithm
[in]infoctchangefstr_contactChange

Definition at line 61 of file fstr_contact.f90.

Here is the caller graph for this function:

◆ fstr_is_matrixstructure_changed()

logical function mcontact::fstr_is_matrixstructure_changed ( type (fstr_info_contactchange)  infoCTChange)
Parameters
infoctchangefstr_contactChange

Definition at line 82 of file fstr_contact.f90.

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

◆ fstr_scan_contact_state()

subroutine mcontact::fstr_scan_contact_state ( integer(kind=kint), intent(in)  cstep,
integer(kind=kint), intent(in)  sub_step,
integer(kind=kint), intent(in)  cont_step,
real(kind=kreal), intent(in)  dt,
integer(kind=kint), intent(in)  ctAlgo,
type( hecmwst_local_mesh ), intent(in)  hecMESH,
type(fstr_solid), intent(inout)  fstrSOLID,
type(fstr_info_contactchange), intent(inout)  infoCTChange,
real(kind=kreal), dimension(:), optional  B 
)

Scanning contact state.

Parameters
[in]cstepcurrent step number
[in]sub_stepcurrent sub-step number
[in]cont_stepcurrent contact step number
[in]ctalgocontact analysis algorithm
[in]hecmeshtype mesh
[in,out]fstrsolidtype fstr_solid
bnodal force residual

Definition at line 212 of file fstr_contact.f90.

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

◆ fstr_scan_contact_state_exp()

subroutine mcontact::fstr_scan_contact_state_exp ( integer(kind=kint), intent(in)  cstep,
type( hecmwst_local_mesh ), intent(in)  hecMESH,
type(fstr_solid), intent(inout)  fstrSOLID,
type(fstr_info_contactchange), intent(inout)  infoCTChange 
)

Scanning contact state.

Parameters
[in]cstepcurrent step number
[in]hecmeshtype mesh
[in,out]fstrsolidtype fstr_solid

Definition at line 350 of file fstr_contact.f90.

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

◆ fstr_set_contact_active()

subroutine mcontact::fstr_set_contact_active ( logical, intent(in)  a)

Definition at line 56 of file fstr_contact.f90.

Here is the caller graph for this function:

◆ fstr_set_contact_penalty()

subroutine mcontact::fstr_set_contact_penalty ( real(kind=kreal), intent(in)  maxv)

Definition at line 45 of file fstr_contact.f90.

Here is the caller graph for this function:

◆ fstr_setup_parancon_contactvalue()

subroutine mcontact::fstr_setup_parancon_contactvalue ( type(hecmwst_local_mesh), intent(in)  hecMESH,
integer(kind=kint), intent(in)  ndof,
real(kind=kreal), dimension(:), intent(inout), pointer  vec,
integer(kind=kint), intent(in)  vtype 
)

Definition at line 764 of file fstr_contact.f90.

Here is the caller graph for this function:

◆ fstr_update_contact0()

subroutine mcontact::fstr_update_contact0 ( type( hecmwst_local_mesh ), intent(in)  hecMESH,
type(fstr_solid), intent(inout)  fstrSOLID,
real(kind=kreal), dimension(:), intent(inout)  B 
)

Update lagrangian multiplier.

Parameters
[in]hecmeshtype mesh
[in,out]fstrsolidtype fstr_solid
[in,out]bnodal force residual

Definition at line 393 of file fstr_contact.f90.

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

◆ fstr_update_contact_multiplier()

subroutine mcontact::fstr_update_contact_multiplier ( type( hecmwst_local_mesh ), intent(in)  hecMESH,
type(fstr_solid), intent(inout)  fstrSOLID,
logical, intent(out)  ctchanged 
)

Update lagrangian multiplier.

Definition at line 418 of file fstr_contact.f90.

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

◆ fstr_update_contact_tangentforce()

subroutine mcontact::fstr_update_contact_tangentforce ( type(fstr_solid), intent(inout)  fstrSOLID)

Update tangent force.

Definition at line 469 of file fstr_contact.f90.

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

◆ fstr_write_mpc()

subroutine mcontact::fstr_write_mpc ( integer(kind=kint), intent(in)  file,
type( hecmwst_mpc ), intent(in)  mpcs 
)

Print out mpc conditions.

Parameters
[in]filefile number
[in]mpcsmpcs to be printed

Definition at line 194 of file fstr_contact.f90.

◆ initialize_contact_output_vectors()

subroutine mcontact::initialize_contact_output_vectors ( type(fstr_solid fstrSOLID,
type(hecmwst_matrix hecMAT 
)
Parameters
fstrsolidtype fstr_solid
hecmattype hecmwST_matrix

Definition at line 565 of file fstr_contact.f90.

Here is the caller graph for this function:

◆ initialize_embed_vectors()

subroutine mcontact::initialize_embed_vectors ( type(fstr_solid fstrSOLID,
type(hecmwst_matrix hecMAT 
)
Parameters
fstrsolidtype fstr_solid
hecmattype hecmwST_matrix

Definition at line 606 of file fstr_contact.f90.

Here is the caller graph for this function:

◆ print_contatct_pair()

subroutine mcontact::print_contatct_pair ( integer(kind=kint), intent(in)  file,
type( hecmwst_contact_pair ), intent(in)  pair 
)

Write out the contact definition read from mesh file.

Definition at line 33 of file fstr_contact.f90.

◆ remove_duplication_tiedcontact()

subroutine mcontact::remove_duplication_tiedcontact ( integer(kind=kint), intent(in)  cstep,
type( hecmwst_local_mesh ), intent(in)  hecMESH,
type(fstr_solid), intent(inout)  fstrSOLID,
type(fstr_info_contactchange), intent(inout)  infoCTChange 
)

Scanning contact state.

Parameters
[in]cstepcurrent step number
[in]hecmeshtype mesh
[in,out]fstrsolidtype fstr_solid

Definition at line 295 of file fstr_contact.f90.

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

◆ setup_contact_elesurf_for_area()

subroutine mcontact::setup_contact_elesurf_for_area ( integer(kind=kint), intent(in)  cstep,
type( hecmwst_local_mesh ), intent(in)  hecMESH,
type(fstr_solid), intent(inout)  fstrSOLID 
)
Parameters
[in]cstepcurrent step number
[in]hecmeshtype mesh
[in,out]fstrsolidtype fstr_solid

Definition at line 616 of file fstr_contact.f90.

Here is the caller graph for this function:

Variable Documentation

◆ bakgnt

real(kind=kreal), dimension(2), save mcontact::bakgnt

1:current average penetration; 2:current relative tangent displacement!

Definition at line 26 of file fstr_contact.f90.

◆ cdotp

real(kind=kreal), save mcontact::cdotp =1.d3

mu=cdotp*maxval

Definition at line 19 of file fstr_contact.f90.

◆ cgn

real(kind=kreal), save mcontact::cgn =1.d-5

convergent condition of penetration

Definition at line 21 of file fstr_contact.f90.

◆ cgt

real(kind=kreal), save mcontact::cgt =1.d-3

convergent condition of relative tangent disp

Definition at line 22 of file fstr_contact.f90.

◆ gnt

real(kind=kreal), dimension(2), save mcontact::gnt

1:current average penetration; 2:current relative tangent displacement

Definition at line 24 of file fstr_contact.f90.

◆ mu

real(kind=kreal), save mcontact::mu =1.d10

penalty, default value

Definition at line 17 of file fstr_contact.f90.

◆ mut

real(kind=kreal), save mcontact::mut =1.d6

penalty along tangent direction

Definition at line 18 of file fstr_contact.f90.

◆ n_contact_mpc

integer(kind=kint), save mcontact::n_contact_mpc

Definition at line 14 of file fstr_contact.f90.