![]() |
FrontISTR
5.9.0
Large-scale structural analysis program with finit element method
|
This module provides geometric calculations for contact. More...
Functions/Subroutines | |
| subroutine | project_point2element (xyz, etype, nn, elemt, reflen, cstate, isin, distclr, ctpos, localclr) |
| This subroutine find the projection of a slave point onto master surface. More... | |
| subroutine | project_point2solidelement (xyz, etype, nn, elemt, reflen, cstate, isin, distclr, ctpos, localclr) |
| This subroutine find the projection of a slave point onto master surface. More... | |
| subroutine | project_point2surfelement (xyz, surf, currpos, cstate, isin, distclr, ctpos, localclr) |
| Wrapper for project_Point2Element that takes tSurfElement structure This subroutine handles element coordinate extraction from tSurfElement. More... | |
| subroutine | getmetrictensor (pos, etype, ele, tensor) |
| This subroutine calculate the metric tensor of a elemental surface. More... | |
| subroutine | dispincrematrix (pos, etype, nnode, ele, tangent, tensor, matrix) |
| This subroutine calculate the relation between global disp and displacement along natural coordinate of master surface supposing penetration is small. More... | |
| subroutine | update_direction (nslave, contact, currpos) |
| This subroutine tracks down next contact position after a finite slide. More... | |
| subroutine, public | cal_node_normal (csurf, isin, surf, currpos, lpos, normal) |
| Calculate averaged nodal normal. More... | |
This module provides geometric calculations for contact.
| subroutine, public m_fstr_contact_geom::cal_node_normal | ( | integer, intent(in) | csurf, |
| integer, intent(in) | isin, | ||
| type(tsurfelement), dimension(:), intent(in) | surf, | ||
| real(kind=kreal), dimension(:), intent(in) | currpos, | ||
| real(kind=kreal), dimension(:), intent(in) | lpos, | ||
| real(kind=kreal), dimension(3), intent(out) | normal | ||
| ) |
Calculate averaged nodal normal.
| [in] | csurf | current surface element |
| [in] | isin | return value from isInsideElement() |
| [in] | surf | surface elements |
| [in] | currpos | current coordinate of each nodes |
| [in] | lpos | local coordinate of contact position |
| [out] | normal | averaged node nomral |
Definition at line 346 of file fstr_contact_geom.f90.
| subroutine m_fstr_contact_geom::dispincrematrix | ( | real(kind=kreal), dimension(2), intent(in) | pos, |
| integer, intent(in) | etype, | ||
| integer, intent(in) | nnode, | ||
| real(kind=kreal), dimension(3,nnode), intent(in) | ele, | ||
| real(kind=kreal), dimension(3,2), intent(out) | tangent, | ||
| real(kind=kreal), dimension(2,2), intent(out) | tensor, | ||
| real(kind=kreal), dimension(2,nnode*3+3), intent(out) | matrix | ||
| ) |
This subroutine calculate the relation between global disp and displacement along natural coordinate of master surface supposing penetration is small.
| [in] | pos | current position(local coordinate) |
| [in] | etype | surface element type |
| [in] | nnode | number of nodes in surface |
| [in] | ele | elemental coordinates |
| [out] | tangent | tangent basis |
| [out] | tensor | metric tensor |
| [out] | matrix | relation between local and global disp increment |
Definition at line 268 of file fstr_contact_geom.f90.
| subroutine m_fstr_contact_geom::getmetrictensor | ( | real(kind=kreal), dimension(2), intent(in) | pos, |
| integer, intent(in) | etype, | ||
| real(kind=kreal), dimension(:,:), intent(in) | ele, | ||
| real(kind=kreal), dimension(2,2), intent(out) | tensor | ||
| ) |
This subroutine calculate the metric tensor of a elemental surface.
| [in] | pos | current position(local coordinate) |
| [in] | etype | surface element type |
| [in] | ele | elemental coordinates |
| [out] | tensor | metric tensor |
Definition at line 250 of file fstr_contact_geom.f90.
| subroutine m_fstr_contact_geom::project_point2element | ( | real(kind=kreal), dimension(3), intent(in) | xyz, |
| integer, intent(in) | etype, | ||
| integer, intent(in) | nn, | ||
| real(kind=kreal), dimension(:,:), intent(in) | elemt, | ||
| real(kind=kreal), intent(in) | reflen, | ||
| type(tcontactstate), intent(inout) | cstate, | ||
| logical, intent(out) | isin, | ||
| real(kind=kreal), intent(in) | distclr, | ||
| real(kind=kreal), dimension(2), optional | ctpos, | ||
| real(kind=kreal), optional | localclr | ||
| ) |
This subroutine find the projection of a slave point onto master surface.
| [in] | xyz | Coordinates of a spacial point, whose projecting point is to be computed |
| [in] | etype | surface element type |
| [in] | nn | number of elemental nodes |
| [in] | elemt | nodes coordinates of surface element |
| [in] | reflen | reference length of surface element |
| [in,out] | cstate | Recorde of contact information |
| [out] | isin | in contact or not |
| [in] | distclr | clearance of contact distance |
| ctpos | curr contact position( natural coord ) | |
| localclr | clearance of contact local coord |
Definition at line 18 of file fstr_contact_geom.f90.
| subroutine m_fstr_contact_geom::project_point2solidelement | ( | real(kind=kreal), dimension(3), intent(in) | xyz, |
| integer, intent(in) | etype, | ||
| integer, intent(in) | nn, | ||
| real(kind=kreal), dimension(:,:), intent(in) | elemt, | ||
| real(kind=kreal), intent(in) | reflen, | ||
| type(tcontactstate), intent(inout) | cstate, | ||
| logical, intent(out) | isin, | ||
| real(kind=kreal), intent(in) | distclr, | ||
| real(kind=kreal), dimension(3), optional | ctpos, | ||
| real(kind=kreal), optional | localclr | ||
| ) |
This subroutine find the projection of a slave point onto master surface.
| [in] | xyz | Coordinates of a spacial point, whose projecting point is to be computed |
| [in] | etype | surface element type |
| [in] | nn | number of elemental nodes |
| [in] | elemt | nodes coordinates of surface element |
| [in] | reflen | reference length of surface element |
| [in,out] | cstate | Recorde of contact information |
| [out] | isin | in contact or not |
| [in] | distclr | clearance of contact distance |
| ctpos | curr contact position( natural coord ) | |
| localclr | clearance of contact local coord |
Definition at line 138 of file fstr_contact_geom.f90.
| subroutine m_fstr_contact_geom::project_point2surfelement | ( | real(kind=kreal), dimension(3), intent(in) | xyz, |
| type(tsurfelement), intent(in) | surf, | ||
| real(kind=kreal), dimension(:), intent(in) | currpos, | ||
| type(tcontactstate), intent(inout) | cstate, | ||
| logical, intent(out) | isin, | ||
| real(kind=kreal), intent(in) | distclr, | ||
| real(kind=kreal), dimension(2), intent(in), optional | ctpos, | ||
| real(kind=kreal), intent(in), optional | localclr | ||
| ) |
Wrapper for project_Point2Element that takes tSurfElement structure This subroutine handles element coordinate extraction from tSurfElement.
| [in] | xyz | coordinates of slave point |
| [in] | surf | surface element structure |
| [in] | currpos | current coordinate of all nodes |
| [in,out] | cstate | contact state |
| [out] | isin | in contact or not |
| [in] | distclr | clearance of contact distance |
| [in] | ctpos | current contact position (natural coord) |
| [in] | localclr | clearance of contact local coord |
Definition at line 224 of file fstr_contact_geom.f90.
| subroutine m_fstr_contact_geom::update_direction | ( | integer, intent(in) | nslave, |
| type( tcontact ), intent(inout) | contact, | ||
| real(kind=kreal), dimension(:), intent(in) | currpos | ||
| ) |
This subroutine tracks down next contact position after a finite slide.
| [in] | nslave | slave node |
| [in,out] | contact | contact info |
| [in] | currpos | current coordinate of each nodes |
checking the contact element of last step
Definition at line 309 of file fstr_contact_geom.f90.