![]() |
FrontISTR
5.9.0
Large-scale structural analysis program with finit element method
|
Contact damping module for CONTACTNEAR state. More...
Functions/Subroutines | |
| pure logical function, public | is_damping_enabled (contact) |
| Check if damping is enabled for a contact pair. More... | |
| pure real(kind=kreal) function, public | contact_damping_weight (g, gact) |
| Compute smoothstep weight w(g) for contact damping. More... | |
| subroutine, public | getdampingstiffness (cstate, surf, alpha, gact, stiff, smoothing_type) |
| Compute damping stiffness matrix K_damp = alpha * w(g) * Bn (x) Bn. More... | |
| subroutine, public | getdampingnodalforce (cstate, surf, ndDu, alpha, gact, ctNForce, ctTForce, smoothing_type) |
| Compute damping nodal force vector f_damp = alpha * w(g) * (Bn^T * ddu) * Bn. More... | |
Contact damping module for CONTACTNEAR state.
Provides smoothstep weight function w(g), and routines to compute damping nodal force and tangent stiffness for NEAR-state slave nodes.
Damping residual: r_damp = alpha * w(g) * (Bn^T * ddu) * Bn Damping stiffness: K_damp = alpha * w(g) * Bn (x) Bn
where g = cstatedistance (frozen during Newton iteration), alpha = DAMP_ALPHA * refStiff, and w(g) is a smoothstep that transitions from 1 at g=0 to 0 at g=gact.
| pure real(kind=kreal) function, public m_fstr_contact_damping::contact_damping_weight | ( | real(kind=kreal), intent(in) | g, |
| real(kind=kreal), intent(in) | gact | ||
| ) |
Compute smoothstep weight w(g) for contact damping.
w(g) = 1 if g <= 0 (penetrating) w(g) = 1 - 3t^2 + 2t^3 if 0 < g < gact (near zone) w(g) = 0 if g >= gact (far away)
where t = g / gact
| [in] | g | signed distance (+ = gap, - = penetration) |
| [in] | gact | activation distance (NEAR_DIST) |
Definition at line 44 of file fstr_contact_damping.f90.
| subroutine, public m_fstr_contact_damping::getdampingnodalforce | ( | type(tcontactstate), intent(in) | cstate, |
| type(tsurfelement), intent(in) | surf, | ||
| real(kind=kreal), dimension(:), intent(in) | ndDu, | ||
| real(kind=kreal), intent(in) | alpha, | ||
| real(kind=kreal), intent(in) | gact, | ||
| real(kind=kreal), dimension(:), intent(out) | ctNForce, | ||
| real(kind=kreal), dimension(:), intent(out) | ctTForce, | ||
| integer(kind=kint), intent(in), optional | smoothing_type | ||
| ) |
Compute damping nodal force vector f_damp = alpha * w(g) * (Bn^T * ddu) * Bn.
Same interface pattern as getContactNodalForce_Alag: compute force vector, caller assembles via assemble_contact_force_residual. ctNForce receives the damping force; ctTForce is zeroed.
| [in] | cstate | contact state (NEAR) |
| [in] | surf | master surface element |
| [in] | nddu | nodal displacement increment (ddisp for this contact element) |
| [in] | alpha | damping coefficient (DAMP_ALPHA * refStiff) |
| [in] | gact | activation distance (NEAR_DIST) |
| [out] | ctnforce | damping force vector (normal-like) |
| [out] | cttforce | tangential force vector (zeroed) |
Definition at line 105 of file fstr_contact_damping.f90.
| subroutine, public m_fstr_contact_damping::getdampingstiffness | ( | type(tcontactstate), intent(in) | cstate, |
| type(tsurfelement), intent(in) | surf, | ||
| real(kind=kreal), intent(in) | alpha, | ||
| real(kind=kreal), intent(in) | gact, | ||
| real(kind=kreal), dimension(:,:), intent(out) | stiff, | ||
| integer(kind=kint), intent(in), optional | smoothing_type | ||
| ) |
Compute damping stiffness matrix K_damp = alpha * w(g) * Bn (x) Bn.
Same interface pattern as getContactStiffness_Alag: compute stiffness matrix, caller assembles via hecmw_mat_ass_elem.
| [in] | cstate | contact state (NEAR) |
| [in] | surf | master surface element |
| [in] | alpha | damping coefficient (DAMP_ALPHA * refStiff) |
| [in] | gact | activation distance (NEAR_DIST) |
| [out] | stiff | damping stiffness matrix |
Definition at line 63 of file fstr_contact_damping.f90.
| pure logical function, public m_fstr_contact_damping::is_damping_enabled | ( | type(tcontact), intent(in) | contact | ) |
Check if damping is enabled for a contact pair.
Definition at line 32 of file fstr_contact_damping.f90.