![]() |
FrontISTR
5.9.0
Large-scale structural analysis program with finit element method
|
Contact processing at assembly level (all pairs in one tContact object) More...
Functions/Subroutines | |
| subroutine, public | calc_contact_pair_refstiff (contact, diag, ndof, hecMESH) |
| Calculate reference stiffness for one contact pair. More... | |
| subroutine | assemble_contact_force_residual (nnode, ndLocal, id_lagrange, ctNForce, ctTForce, conMAT) |
| Assemble contact nodal force into residual vector (conMATB). More... | |
| subroutine | assemble_contact_force_output (nnode, ndLocal, ctNForce, ctTForce, cont_nforce, cont_fric) |
| Accumulate contact nodal force into output arrays (CONT_NFORCE/CONT_FRIC). More... | |
| subroutine | update_contact_multiplier (ctAlgo, contact, coord, disp, ddisp, fcoeff, hecMESH, hecLagMAT, gnt, ctchanged) |
| This subroutine update lagrangian multiplier and the distance between contacting nodes. More... | |
| subroutine | update_tied_multiplier (contact, disp, ddisp, ctchanged) |
| This subroutine update lagrangian multiplier and the distance between contacting nodes. More... | |
| subroutine | update_contact_tangentforce (contact) |
| subroutine | calcu_contact_stiffness_nodesurf (ctAlgo, contact, coord, disp, iter, lagrange_array, conMAT, hecLagMAT) |
| This subroutine calculates contact stiffness for each contact pair and assembles it into global stiffness matrix. More... | |
| subroutine | calcu_contact_ndforce_nodesurf (purpose, ctAlgo, contact, coord, disp, ddisp, lagrange_array, conMAT, CONT_NFORCE, CONT_FRIC, hecLagMAT) |
| This subroutine calculates contact nodal force for each contact pair and assembles it into contact matrix and/or force arrays. When purpose == kctForResidual, forces are assembled into conMATB. When purpose == kctForOutput, forces are stored in CONT_NFORCE/CONT_FRIC using multiplier only (no penalty). More... | |
Contact processing at assembly level (all pairs in one tContact object)
| subroutine m_fstr_contact_assembly::assemble_contact_force_output | ( | integer(kind=kint), intent(in) | nnode, |
| integer(kind=kint), dimension(nnode + 1), intent(in) | ndLocal, | ||
| real(kind=kreal), dimension((nnode+1)*3+1), intent(in) | ctNForce, | ||
| real(kind=kreal), dimension((nnode+1)*3+1), intent(in) | ctTForce, | ||
| real(kind=kreal), dimension(:), intent(inout), pointer | cont_nforce, | ||
| real(kind=kreal), dimension(:), intent(inout), optional, pointer | cont_fric | ||
| ) |
Accumulate contact nodal force into output arrays (CONT_NFORCE/CONT_FRIC).
| [in] | nnode | number of master nodes |
| [in] | ndlocal | global node numbers |
| [in] | ctnforce | normal contact force vector |
| [in] | cttforce | tangential contact force vector |
| [in,out] | cont_nforce | output normal force |
| [in,out] | cont_fric | output friction force |
Definition at line 92 of file fstr_contact_assembly.f90.
| subroutine m_fstr_contact_assembly::assemble_contact_force_residual | ( | integer(kind=kint), intent(in) | nnode, |
| integer(kind=kint), dimension(nnode + 1), intent(in) | ndLocal, | ||
| integer(kind=kint), intent(in) | id_lagrange, | ||
| real(kind=kreal), dimension((nnode+1)*3+1), intent(in) | ctNForce, | ||
| real(kind=kreal), dimension((nnode+1)*3+1), intent(in) | ctTForce, | ||
| type(hecmwst_matrix), intent(inout) | conMAT | ||
| ) |
Assemble contact nodal force into residual vector (conMATB).
| [in] | nnode | number of master nodes |
| [in] | ndlocal | global node numbers |
| [in] | id_lagrange | Lagrange multiplier index (0 if none) |
| [in] | ctnforce | normal contact force vector |
| [in] | cttforce | tangential contact force vector |
| [in,out] | conmat | contact matrix |
Definition at line 69 of file fstr_contact_assembly.f90.
| subroutine, public m_fstr_contact_assembly::calc_contact_pair_refstiff | ( | type(tcontact), intent(inout) | contact, |
| real(kind=kreal), dimension(:), intent(in) | diag, | ||
| integer(kind=kint), intent(in) | ndof, | ||
| type(hecmwst_local_mesh), intent(in) | hecMESH | ||
| ) |
Calculate reference stiffness for one contact pair.
| [in,out] | contact | contact pair |
| [in] | diag | diagonal vector (size = ndof * np) |
| [in] | ndof | degrees of freedom |
| [in] | hecmesh | mesh |
Definition at line 21 of file fstr_contact_assembly.f90.
| subroutine m_fstr_contact_assembly::calcu_contact_ndforce_nodesurf | ( | integer(kind=kint), intent(in) | purpose, |
| integer(kind=kint), intent(in) | ctAlgo, | ||
| type( tcontact ), intent(inout) | contact, | ||
| real(kind=kreal), dimension(:), intent(in) | coord, | ||
| real(kind=kreal), dimension(:), intent(in) | disp, | ||
| real(kind=kreal), dimension(:), intent(in) | ddisp, | ||
| real(kind=kreal), dimension(:), intent(in) | lagrange_array, | ||
| type(hecmwst_matrix), intent(inout) | conMAT, | ||
| real(kind=kreal), dimension(:), pointer | CONT_NFORCE, | ||
| real(kind=kreal), dimension(:), pointer | CONT_FRIC, | ||
| type(hecmwst_matrix_lagrange), intent(in) | hecLagMAT | ||
| ) |
This subroutine calculates contact nodal force for each contact pair and assembles it into contact matrix and/or force arrays. When purpose == kctForResidual, forces are assembled into conMATB. When purpose == kctForOutput, forces are stored in CONT_NFORCE/CONT_FRIC using multiplier only (no penalty).
| [in] | purpose | kctForResidual or kctForOutput |
| [in] | ctalgo | contact analysis algorithm |
| [in,out] | contact | contact info |
| [in] | coord | mesh coordinate |
| [in] | disp | disp till current step |
| [in] | ddisp | disp till current substep |
| [in] | lagrange_array | Lagrange multiplier array |
| [in,out] | conmat | contact matrix |
| cont_nforce | contact normal force | |
| cont_fric | contact friction force | |
| [in] | heclagmat | Lagrange matrix |
Definition at line 342 of file fstr_contact_assembly.f90.
| subroutine m_fstr_contact_assembly::calcu_contact_stiffness_nodesurf | ( | integer(kind=kint), intent(in) | ctAlgo, |
| type(tcontact), intent(inout) | contact, | ||
| real(kind=kreal), dimension(:), intent(in) | coord, | ||
| real(kind=kreal), dimension(:), intent(in) | disp, | ||
| integer(kind=kint), intent(in) | iter, | ||
| real(kind=kreal), dimension(:), intent(in) | lagrange_array, | ||
| type(hecmwst_matrix), intent(inout) | conMAT, | ||
| type(hecmwst_matrix_lagrange), intent(inout) | hecLagMAT | ||
| ) |
This subroutine calculates contact stiffness for each contact pair and assembles it into global stiffness matrix.
| [in] | ctalgo | contact analysis algorithm |
| [in,out] | contact | contact info |
| [in] | coord | mesh coordinate |
| [in] | disp | displacement |
| [in] | iter | iteration number |
| [in] | lagrange_array | Lagrange multiplier array |
| [in,out] | conmat | contact stiffness matrix |
| [in,out] | heclagmat | Lagrange matrix |
Definition at line 252 of file fstr_contact_assembly.f90.
| subroutine m_fstr_contact_assembly::update_contact_multiplier | ( | integer(kind=kint), intent(in) | ctAlgo, |
| type( tcontact ), intent(inout) | contact, | ||
| real(kind=kreal), dimension(:), intent(in) | coord, | ||
| real(kind=kreal), dimension(:), intent(in) | disp, | ||
| real(kind=kreal), dimension(:), intent(in) | ddisp, | ||
| real(kind=kreal), intent(in) | fcoeff, | ||
| type(hecmwst_local_mesh), intent(in) | hecMESH, | ||
| type(hecmwst_matrix_lagrange), intent(in) | hecLagMAT, | ||
| real(kind=kreal), dimension(2), intent(out) | gnt, | ||
| logical, intent(inout) | ctchanged | ||
| ) |
This subroutine update lagrangian multiplier and the distance between contacting nodes.
| [in] | ctalgo | contact algorithm |
| [in,out] | contact | contact info |
| [in] | coord | mesh coordinate |
| [in] | disp | disp till current step |
| [in] | ddisp | disp till current substep |
| [in] | fcoeff | frictional coeff |
| [in] | hecmesh | mesh for allreduce |
| [in] | heclagmat | Lagrange matrix |
| [out] | gnt | convergency information |
| [in,out] | ctchanged | if contact state changes |
Definition at line 113 of file fstr_contact_assembly.f90.
| subroutine m_fstr_contact_assembly::update_contact_tangentforce | ( | type( tcontact ), intent(inout) | contact | ) |
| [in,out] | contact | contact info |
Definition at line 233 of file fstr_contact_assembly.f90.
| subroutine m_fstr_contact_assembly::update_tied_multiplier | ( | type( tcontact ), intent(inout) | contact, |
| real(kind=kreal), dimension(:), intent(in) | disp, | ||
| real(kind=kreal), dimension(:), intent(in) | ddisp, | ||
| logical, intent(inout) | ctchanged | ||
| ) |
This subroutine update lagrangian multiplier and the distance between contacting nodes.
| [in,out] | contact | contact info |
| [in] | disp | disp till current step |
| [in] | ddisp | disp till current substep |
| [in,out] | ctchanged | if contact state changes |
Definition at line 182 of file fstr_contact_assembly.f90.