![]() |
FrontISTR
5.9.0
Large-scale structural analysis program with finit element method
|
Common utilities for contact element calculations. More...
Functions/Subroutines | |
| subroutine, public | computetm_tt (ctState, tSurf, fcoeff, Tm, Tt) |
| Compute Tm (relative displacement mapping) and optionally Tt (tangential mapping) This subroutine constructs the mapping matrices based on contact_disc.md section 10: Tm = [I_3; -N_1*I_3; -N_2*I_3; ...; -N_n*I_3] (size: 3 x 3*(nnode+1)) Tt = Pt * Tm where Pt = I_3 - n x n (computed only if fcoeff /= 0) More... | |
| subroutine, public | computecontactmaps_alag (ctState, tSurf, ele, Bn, metric, Ht, Gt) |
| Compute contact maps for ALag method (normal distribution and tangential displacement maps) This subroutine constructs the mapping matrices for ALag contact using computeTm_Tt as the core: Bn: normal force distribution vector [n; -N_1*n; -N_2*n; ...; -N_n*n] (size: nnode*3+3) Ht: tangent displacement map from DispIncreMatrix (size: 2 x nnode*3+3) Gt: covariant displacement map Gt = metric * Ht (size: 2 x nnode*3+3) Where: dxy = Gt * edisp, f3 = Ht^T * fric. More... | |
| subroutine, public | computefrictionforce_alag (ctState, fcoeff, lambda_n, metric, Ht, Gt, edisp, edof, ctTForce, mut, update_multiplier, slave_id, ctchanged, norm_trial, alpha, that, jump_ratio) |
| Compute friction force for ALag method Given tangent maps Ht, Gt and displacement increment, compute the friction force vector. Handles both stick and slip states. More... | |
| subroutine, public | gettrialfricforceandcheckfricstate (ctstate, fcoeff, tPenalty, lagrange) |
| This subroutine calculates trial friction force and checks friction state. More... | |
| subroutine, public | update_tangentforce (etype, nn, elemt0, elemt, cstate) |
| This subroutine find the projection of a slave point onto master surface. More... | |
Common utilities for contact element calculations.
| subroutine, public m_fstr_contact_elem_common::computecontactmaps_alag | ( | type(tcontactstate), intent(in) | ctState, |
| type(tsurfelement), intent(in) | tSurf, | ||
| real(kind=kreal), dimension(:,:), intent(in) | ele, | ||
| real(kind=kreal), dimension(:), intent(out) | Bn, | ||
| real(kind=kreal), dimension(2,2), intent(out) | metric, | ||
| real(kind=kreal), dimension(:,:), intent(out) | Ht, | ||
| real(kind=kreal), dimension(:,:), intent(out) | Gt | ||
| ) |
Compute contact maps for ALag method (normal distribution and tangential displacement maps) This subroutine constructs the mapping matrices for ALag contact using computeTm_Tt as the core: Bn: normal force distribution vector [n; -N_1*n; -N_2*n; ...; -N_n*n] (size: nnode*3+3) Ht: tangent displacement map from DispIncreMatrix (size: 2 x nnode*3+3) Gt: covariant displacement map Gt = metric * Ht (size: 2 x nnode*3+3) Where: dxy = Gt * edisp, f3 = Ht^T * fric.
This routine now uses computeTm_Tt to generate Bn consistently, enabling future smoothing extensions in one place.
| [in] | ctstate | contact state (contains lpos and direction) |
| [in] | tsurf | surface element structure |
| [in] | ele | master node coords (3, nnode): coord+disp |
| [out] | bn | normal distribution vector (nnode*3+3) |
| [out] | metric | metric tensor |
| [out] | ht | tangent displacement map (2, nnode*3+3) |
| [out] | gt | covariant displacement map (2, nnode*3+3) |
Definition at line 95 of file fstr_contact_elem_common.f90.
| subroutine, public m_fstr_contact_elem_common::computefrictionforce_alag | ( | type(tcontactstate), intent(inout) | ctState, |
| real(kind=kreal), intent(in) | fcoeff, | ||
| real(kind=kreal), intent(in) | lambda_n, | ||
| real(kind=kreal), dimension(2,2), intent(in) | metric, | ||
| real(kind=kreal), dimension(:,:), intent(in) | Ht, | ||
| real(kind=kreal), dimension(:,:), intent(in) | Gt, | ||
| real(kind=kreal), dimension(:), intent(in) | edisp, | ||
| integer(kind=kint), intent(in) | edof, | ||
| real(kind=kreal), dimension(:), intent(out) | ctTForce, | ||
| real(kind=kreal), intent(in) | mut, | ||
| logical, intent(in), optional | update_multiplier, | ||
| integer(kind=kint), intent(in), optional | slave_id, | ||
| logical, intent(inout), optional | ctchanged, | ||
| real(kind=kreal), intent(out), optional | norm_trial, | ||
| real(kind=kreal), intent(out), optional | alpha, | ||
| real(kind=kreal), dimension(2), intent(out), optional | that, | ||
| real(kind=kreal), intent(out), optional | jump_ratio | ||
| ) |
Compute friction force for ALag method Given tangent maps Ht, Gt and displacement increment, compute the friction force vector. Handles both stick and slip states.
Slip criterion uses 2D local norm: ||fric||_g = sqrt(fric^T * inv(metric) * fric) This ensures consistency between force calculation and stiffness matrix.
Input: ctState: contact state (contains multiplier(2:3) for tangent Lagrange multiplier) fcoeff: friction coefficient lambda_n: normal Lagrange multiplier (for slip scaling) metric: metric tensor (2x2) for 2D local tangent space Ht: tangent displacement map (2 x edof) Gt: covariant displacement map (2 x edof) edisp: displacement increment vector (edof) edof: element DOF size update_multiplier: (optional) if true, update ctStatemultiplier(2:3) and check state change slave_id: (optional) slave node ID for debug print
Output: ctTForce: friction force vector (edof) ctState: (if update_multiplier=true) updated multiplier(2:3) and state ctchanged: (optional) flag indicating if contact state changed norm_trial: (optional) norm of fric_trial in metric space alpha: (optional) projection ratio = min(1, r/norm_trial), r=fcoeff*lambda_n that: (optional) normalized direction = fric_trial/norm_trial (2D)
| [in,out] | ctstate | contact state (inout for multiplier update) |
| [in] | fcoeff | friction coefficient |
| [in] | lambda_n | normal Lagrange multiplier |
| [in] | metric | metric tensor (2x2) |
| [in] | ht | tangent displacement map (2, edof) |
| [in] | gt | covariant displacement map (2, edof) |
| [in] | edisp | displacement increment |
| [in] | edof | element DOF size |
| [out] | cttforce | friction force vector |
| [in] | mut | tangential penalty parameter |
| [in] | update_multiplier | if true, update multiplier and check state |
| [in] | slave_id | slave node ID for debug print |
| [in,out] | ctchanged | flag for state change |
| [out] | norm_trial | norm of fric_trial (for stiffness) |
| [out] | alpha | projection ratio (for stiffness) |
| [out] | that | normalized direction (for stiffness) |
| [out] | jump_ratio | stick trial / slip limit force ratio |
Definition at line 175 of file fstr_contact_elem_common.f90.
| subroutine, public m_fstr_contact_elem_common::computetm_tt | ( | type(tcontactstate), intent(in) | ctState, |
| type(tsurfelement), intent(in) | tSurf, | ||
| real(kind=kreal), intent(in) | fcoeff, | ||
| real(kind=kreal), dimension(3, 3*(l_max_surface_node+1)), intent(out) | Tm, | ||
| real(kind=kreal), dimension(3, 3*(l_max_surface_node+1)), intent(out) | Tt | ||
| ) |
Compute Tm (relative displacement mapping) and optionally Tt (tangential mapping) This subroutine constructs the mapping matrices based on contact_disc.md section 10: Tm = [I_3; -N_1*I_3; -N_2*I_3; ...; -N_n*I_3] (size: 3 x 3*(nnode+1)) Tt = Pt * Tm where Pt = I_3 - n x n (computed only if fcoeff /= 0)
| [in] | ctstate | contact state (contains lpos and direction) |
| [in] | tsurf | surface element structure |
| [in] | fcoeff | friction coefficient (if 0, Tt not computed) |
| [out] | tm | relative displacement mapping matrix |
| [out] | tt | tangential mapping matrix |
Definition at line 25 of file fstr_contact_elem_common.f90.
| subroutine, public m_fstr_contact_elem_common::gettrialfricforceandcheckfricstate | ( | type(tcontactstate) | ctstate, |
| real(kind=kreal) | fcoeff, | ||
| real(kind=kreal) | tPenalty, | ||
| real(kind=kreal) | lagrange | ||
| ) |
This subroutine calculates trial friction force and checks friction state.
| ctstate | type tContactState |
| tpenalty | friction coefficient; tangential penalty |
| lagrange | value of Lagrange multiplier |
Definition at line 293 of file fstr_contact_elem_common.f90.
| subroutine, public m_fstr_contact_elem_common::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.
| [in] | etype | surface element type |
| [in] | nn | number of elemental nodes |
| [in] | elemt0 | nodes coordinates of surface element at t |
| [in] | elemt | nodes coordinates of surface element at t+dt |
| [in,out] | cstate | Recorde of contact information |
Definition at line 317 of file fstr_contact_elem_common.f90.