![]() |
FrontISTR
5.9.0
Large-scale structural analysis program with finit element method
|
Top-level contact analysis module (System level) More...
Functions/Subroutines | |
| subroutine | fstr_addcontactstiffness (cstep, ctAlgo, iter, hecMESH, conMAT, hecLagMAT, fstrSOLID) |
| subroutine | fstr_update_ndforce_contact (cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT) |
| Compute contact forces for residual vector (conMATB). More... | |
| subroutine | fstr_calc_contact_output_force (cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT) |
| Compute contact forces for output (CONT_NFORCE/CONT_FRIC). More... | |
| subroutine | fstr_contact_ndforce_core (purpose, cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT) |
| Core routine: compute contact nodal forces for all contact/embed pairs. purpose == kctForResidual: assemble into conMATB purpose == kctForOutput: store into CONT_NFORCE/CONT_FRIC (zero-cleared first) More... | |
| subroutine | fstr_calc_contact_refstiff (cstep, hecMESH, hecMAT, fstrSOLID) |
| Calculate reference stiffness for all contact pairs (System Level) More... | |
| subroutine | fstr_scan_contact_state (cstep, sub_step, cont_step, dt, ctAlgo, hecMESH, fstrSOLID, infoCTChange) |
| Scanning contact state. More... | |
| subroutine | fstr_scan_contact_state_exp (cstep, hecMESH, fstrSOLID, infoCTChange) |
| Scanning contact state. More... | |
| 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_update_contact_multiplier (cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, ctchanged) |
| subroutine | fstr_update_contact_tangentforce (cstep, fstrSOLID) |
| Update tangent force. More... | |
| subroutine, public | fstr_update_contact_state_vectors (fstrSOLID, dt) |
| Update contact state output vectors for all contacts. More... | |
Top-level contact analysis module (System level)
Level 1: Element Level (1 slave-master pair) m_fstr_contact_element:
Level 2: Contact Object Level (all pairs in one tContact) m_fstr_contact_assembly:
Level 3: System Level (all contact objects) mContact (this module):
Supporting modules: m_fstr_contact_mpc: MPC (Multi-Point Constraint) processing m_fstr_contact_output: Output vector initialization and processing
| subroutine mcontact::fstr_addcontactstiffness | ( | integer(kind=kint) | cstep, |
| integer(kind=kint) | ctAlgo, | ||
| integer(kind=kint) | iter, | ||
| type(hecmwst_local_mesh) | hecMESH, | ||
| type(hecmwst_matrix) | conMAT, | ||
| type(hecmwst_matrix_lagrange) | hecLagMAT, | ||
| type(fstr_solid) | fstrSOLID | ||
| ) |
| cstep | current loading step |
| ctalgo | contact algorithm type |
| hecmesh | type hecmwST_local_mesh |
| conmat | type hecmwST_matrix |
| fstrsolid | type fstr_solid |
| heclagmat | type hecmwST_matrix_lagrange |
Definition at line 59 of file fstr_contact.f90.
| subroutine mcontact::fstr_calc_contact_output_force | ( | integer(kind=kint), intent(in) | cstep, |
| integer(kind=kint), intent(in) | ctAlgo, | ||
| type(hecmwst_local_mesh) | hecMESH, | ||
| type(hecmwst_matrix_lagrange) | hecLagMAT, | ||
| type(fstr_solid) | fstrSOLID, | ||
| type(hecmwst_matrix) | conMAT | ||
| ) |
Compute contact forces for output (CONT_NFORCE/CONT_FRIC).
Definition at line 109 of file fstr_contact.f90.
| subroutine mcontact::fstr_calc_contact_refstiff | ( | integer(kind=kint), intent(in) | cstep, |
| type(hecmwst_local_mesh) | hecMESH, | ||
| type(hecmwst_matrix) | hecMAT, | ||
| type(fstr_solid) | fstrSOLID | ||
| ) |
Calculate reference stiffness for all contact pairs (System Level)
Definition at line 160 of file fstr_contact.f90.
| subroutine mcontact::fstr_contact_ndforce_core | ( | integer(kind=kint), intent(in) | purpose, |
| integer(kind=kint), intent(in) | cstep, | ||
| integer(kind=kint), intent(in) | ctAlgo, | ||
| type(hecmwst_local_mesh) | hecMESH, | ||
| type(hecmwst_matrix_lagrange) | hecLagMAT, | ||
| type(fstr_solid) | fstrSOLID, | ||
| type(hecmwst_matrix) | conMAT | ||
| ) |
Core routine: compute contact nodal forces for all contact/embed pairs. purpose == kctForResidual: assemble into conMATB purpose == kctForOutput: store into CONT_NFORCE/CONT_FRIC (zero-cleared first)
Definition at line 124 of file fstr_contact.f90.
| logical function mcontact::fstr_is_contact_active |
| 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 | ||
| ) |
| [in] | ctalgo | contact analysis algorithm |
| [in] | infoctchange | fstr_contactChange |
Definition at line 320 of file fstr_contact.f90.
| logical function mcontact::fstr_is_matrixstructure_changed | ( | type (fstr_info_contactchange) | infoCTChange | ) |
| infoctchange | fstr_contactChange |
Definition at line 333 of file fstr_contact.f90.
| 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 | ||
| ) |
Scanning contact state.
| [in] | cstep | current step number |
| [in] | sub_step | current sub-step number |
| [in] | cont_step | current contact step number |
| [in] | ctalgo | contact analysis algorithm |
| [in] | hecmesh | type mesh |
| [in,out] | fstrsolid | type fstr_solid |
Definition at line 195 of file fstr_contact.f90.
| 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.
| [in] | cstep | current step number |
| [in] | hecmesh | type mesh |
| [in,out] | fstrsolid | type fstr_solid |
Definition at line 269 of file fstr_contact.f90.
| subroutine mcontact::fstr_set_contact_active | ( | logical, intent(in) | a | ) |
| subroutine mcontact::fstr_update_contact_multiplier | ( | integer(kind=kint), intent(in) | cstep, |
| integer(kind=kint), intent(in) | ctAlgo, | ||
| type( hecmwst_local_mesh ), intent(in) | hecMESH, | ||
| type(hecmwst_matrix_lagrange), intent(in) | hecLagMAT, | ||
| type(fstr_solid), intent(inout) | fstrSOLID, | ||
| logical, intent(out) | ctchanged | ||
| ) |
Definition at line 340 of file fstr_contact.f90.
| subroutine, public mcontact::fstr_update_contact_state_vectors | ( | type(fstr_solid), intent(inout) | fstrSOLID, |
| real(kind=kreal), intent(in) | dt | ||
| ) |
Update contact state output vectors for all contacts.
| [in] | dt | time increment |
Definition at line 394 of file fstr_contact.f90.
| subroutine mcontact::fstr_update_contact_tangentforce | ( | integer(kind=kint), intent(in) | cstep, |
| type(fstr_solid), intent(inout) | fstrSOLID | ||
| ) |
Update tangent force.
Definition at line 378 of file fstr_contact.f90.
| subroutine mcontact::fstr_update_ndforce_contact | ( | integer(kind=kint), intent(in) | cstep, |
| integer(kind=kint), intent(in) | ctAlgo, | ||
| type(hecmwst_local_mesh) | hecMESH, | ||
| type(hecmwst_matrix_lagrange) | hecLagMAT, | ||
| type(fstr_solid) | fstrSOLID, | ||
| type(hecmwst_matrix) | conMAT | ||
| ) |
Compute contact forces for residual vector (conMATB).
Definition at line 96 of file fstr_contact.f90.