FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_contact_elem_common.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
7  use hecmw
8  use elementinfo
9  use mcontactdef
11  implicit none
12 
13  public :: computetm_tt
14  public :: computecontactmaps_alag
17  public :: update_tangentforce
18 
19 contains
20 
25  subroutine computetm_tt(ctState, tSurf, fcoeff, Tm, Tt)
26  implicit none
27 
28  ! Input arguments
29  type(tcontactstate), intent(in) :: ctstate
30  type(tsurfelement), intent(in) :: tsurf
31  real(kind=kreal), intent(in) :: fcoeff
32 
33  ! Output arguments
34  real(kind=kreal), intent(out) :: tm(3, 3*(l_max_surface_node+1))
35  real(kind=kreal), intent(out) :: tt(3, 3*(l_max_surface_node+1))
36 
37  ! Local variables
38  integer(kind=kint) :: i, j
39  integer(kind=kint) :: nnode
40  real(kind=kreal) :: shapefunc(l_max_surface_node)
41  real(kind=kreal) :: normal(3)
42  real(kind=kreal) :: pt(3,3)
43 
44  tm = 0.0d0
45  tt = 0.0d0
46 
47  nnode = size(tsurf%nodes)
48 
49  call getshapefunc(tsurf%etype, ctstate%lpos(:), shapefunc)
50 
51  ! Get normal vector from contact state
52  normal(1:3) = ctstate%direction(1:3)
53 
54  ! Construct Tm = [I_3; -N_1*I_3; -N_2*I_3; ...; -N_n*I_3]
55  ! First block (slave node): Identity matrix I_3
56  tm(1,1) = 1.0d0
57  tm(2,2) = 1.0d0
58  tm(3,3) = 1.0d0
59 
60  ! Remaining blocks (master nodes): -N_i * I_3
61  do i = 1, nnode
62  tm(1, i*3+1) = -shapefunc(i)
63  tm(2, i*3+2) = -shapefunc(i)
64  tm(3, i*3+3) = -shapefunc(i)
65  enddo
66 
67  ! Compute Tt only if friction coefficient is non-zero
68  if (fcoeff /= 0.0d0) then
69  ! Construct tangential projection operator Pt = I_3 - n x n
70  do i = 1, 3
71  do j = 1, 3
72  if (i == j) then
73  pt(i,j) = 1.0d0 - normal(i)*normal(j)
74  else
75  pt(i,j) = -normal(i)*normal(j)
76  endif
77  enddo
78  enddo
79 
80  ! Compute Tt = Pt * Tm using matrix multiplication
81  tt(1:3, 1:3*(l_max_surface_node+1)) = matmul(pt(1:3,1:3), tm(1:3, 1:3*(l_max_surface_node+1)))
82  endif
83 
84  end subroutine computetm_tt
85 
95  subroutine computecontactmaps_alag(ctState, tSurf, ele, &
96  Bn, metric, Ht, Gt)
97  implicit none
98 
99  ! Input arguments
100  type(tcontactstate), intent(in) :: ctstate
101  type(tsurfelement), intent(in) :: tsurf
102  real(kind=kreal), intent(in) :: ele(:,:)
103 
104  ! Output arguments
105  real(kind=kreal), intent(out) :: bn(:)
106  real(kind=kreal), intent(out) :: metric(2,2)
107  real(kind=kreal), intent(out) :: ht(:,:)
108  real(kind=kreal), intent(out) :: gt(:,:)
109 
110  ! Local variables
111  integer(kind=kint) :: nnode, ndof
112  real(kind=kreal) :: normal(3)
113  real(kind=kreal) :: tangent(3,2)
114  real(kind=kreal) :: tm(3, 3*(l_max_surface_node+1))
115  real(kind=kreal) :: tt(3, 3*(l_max_surface_node+1))
116 
117  ! Initialize arrays to zero for safety
118  bn(:) = 0.0d0
119  metric(:,:) = 0.0d0
120  ht(:,:) = 0.0d0
121  gt(:,:) = 0.0d0
122  tm(:,:) = 0.0d0
123  tt(:,:) = 0.0d0
124 
125  nnode = size(tsurf%nodes)
126  ndof = nnode*3 + 3
127 
128  ! Call computeTm_Tt as the core routine (fcoeff=0 since we don't need Tt for ALag)
129  call computetm_tt(ctstate, tsurf, 0.0d0, tm, tt)
130 
131  ! Get normal vector from contact state
132  normal(1:3) = ctstate%direction(1:3)
133 
134  ! Construct normal distribution vector Bn from Tm
135  ! Bn = Tm^T * normal
136  bn(1:ndof) = matmul(transpose(tm(1:3, 1:ndof)), normal(1:3))
137 
138  ! Call DispIncreMatrix to get tangent displacement map and metric
139  ! DispIncreMatrix returns: tangent basis, metric tensor, and displacement increment matrix
140  call dispincrematrix(ctstate%lpos(1:2), tsurf%etype, nnode, ele(1:3,1:nnode), &
141  tangent, metric, ht(1:2,1:ndof))
142 
143  ! Compute covariant displacement map: Gt = metric * Ht
144  ! This gives: dxy = Gt * edisp = metric * (Ht * edisp)
145  gt(1:2,1:ndof) = matmul(metric(1:2,1:2), ht(1:2,1:ndof))
146 
147  end subroutine computecontactmaps_alag
148 
175  subroutine computefrictionforce_alag(ctState, fcoeff, lambda_n, metric, Ht, Gt, edisp, edof, ctTForce, &
176  mut, update_multiplier, slave_id, ctchanged, &
177  norm_trial, alpha, that, jump_ratio)
178  implicit none
179 
180  type(tcontactstate), intent(inout) :: ctstate
181  real(kind=kreal), intent(in) :: fcoeff
182  real(kind=kreal), intent(in) :: lambda_n
183  real(kind=kreal), intent(in) :: metric(2,2)
184  real(kind=kreal), intent(in) :: ht(:,:)
185  real(kind=kreal), intent(in) :: gt(:,:)
186  real(kind=kreal), intent(in) :: edisp(:)
187  integer(kind=kint), intent(in) :: edof
188  real(kind=kreal), intent(out) :: cttforce(:)
189  real(kind=kreal), intent(in) :: mut
190  logical, optional, intent(in) :: update_multiplier
191  integer(kind=kint), optional, intent(in) :: slave_id
192  logical, optional, intent(inout) :: ctchanged
193  real(kind=kreal), optional, intent(out) :: norm_trial
194  real(kind=kreal), optional, intent(out) :: alpha
195  real(kind=kreal), optional, intent(out) :: that(2)
196  real(kind=kreal), optional, intent(out) :: jump_ratio
197 
198  real(kind=kreal) :: dxy(2)
199  real(kind=kreal) :: fric(2)
200  real(kind=kreal) :: f3(edof)
201  real(kind=kreal) :: invmetric(2,2)
202  real(kind=kreal) :: det
203  real(kind=kreal) :: norm_fric
204  real(kind=kreal) :: tmp(2)
205  real(kind=kreal) :: alpha_local
206  logical :: do_update
207 
208  cttforce = 0.0d0
209  do_update = .false.
210  if(present(jump_ratio)) jump_ratio = 0.0d0
211  if(present(update_multiplier)) do_update = update_multiplier
212 
213  ! Compute inverse of metric tensor (2x2)
214  det = metric(1,1)*metric(2,2) - metric(1,2)*metric(2,1)
215  if( abs(det) < 1.0d-20 ) then
216  ! Degenerate metric: set friction to zero
217  cttforce = 0.0d0
218  return
219  endif
220  invmetric(1,1) = metric(2,2) / det
221  invmetric(2,2) = metric(1,1) / det
222  invmetric(1,2) = -metric(1,2) / det
223  invmetric(2,1) = -metric(2,1) / det
224 
225  ! Compute tangent displacement using Gt map: dxy = Gt * edisp
226  dxy(1:2) = matmul( gt(1:2,1:edof), edisp(1:edof) )
227 
228  ! Trial friction force: fric = λ_t + ρ_t * dxy
229  fric(1:2) = ctstate%multiplier(2:3) + mut*dxy(1:2)
230 
231  ! Compute 2D local norm: ||fric||_g = sqrt(fric^T * inv(metric) * fric)
232  tmp(1:2) = matmul(invmetric(1:2,1:2), fric(1:2))
233  norm_fric = dsqrt( dot_product(fric(1:2), tmp(1:2)) )
234 
235  ! Return projection information for stiffness (optional outputs)
236  if(present(norm_trial)) norm_trial = norm_fric
237  if(present(that) .and. norm_fric > 1.0d-20) then
238  that(1:2) = fric(1:2) / norm_fric
239  else if(present(that)) then
240  that(1:2) = 0.0d0
241  endif
242 
243  ! Compute projection ratio: alpha = min(1, r/norm_trial)
244  if( lambda_n > 0.0d0 .and. norm_fric > 1.0d-20 ) then
245  alpha_local = min(1.0d0, fcoeff*lambda_n / norm_fric)
246  else
247  alpha_local = 0.0d0
248  endif
249  if(present(alpha)) alpha = alpha_local
250 
251  ! Check stick/slip state using projection ratio
252  if( do_update .and. lambda_n > 0.0d0 ) then
253  ! Update mode: check state and update multiplier
254  if( alpha_local < 1.0d0 ) then
255  ! Slip state: projected to friction cone
256  if( ctstate%state == contactstick ) then
257  if(present(ctchanged)) ctchanged = .true.
258  if(present(slave_id)) print *, "Node", slave_id, "to slip state", norm_fric, fcoeff*lambda_n
259  if(present(jump_ratio) .and. lambda_n > 1.0e-20 .and. norm_fric > 1.0e-20) &
260  & jump_ratio = norm_fric / (fcoeff*lambda_n)
261  endif
262  ctstate%state = contactslip
263  fric(1:2) = fric(1:2) * alpha_local
264  else
265  ! Stick state
266  if( ctstate%state == contactslip ) then
267  if(present(ctchanged)) ctchanged = .true.
268  if(present(slave_id)) print *, "Node", slave_id, "to stick state", norm_fric, fcoeff*lambda_n
269  endif
270  ctstate%state = contactstick
271  endif
272  ! Update multiplier
273  ctstate%multiplier(2:3) = fric(1:2)
274  else if( lambda_n <= 0.0d0 ) then
275  ! No normal force: no friction
276  fric(1:2) = 0.0d0
277  else
278  ! Read-only mode: just scale if alpha < 1
279  if( alpha_local < 1.0d0 ) then
280  fric(1:2) = fric(1:2) * alpha_local
281  endif
282  endif
283 
284  ! Compute friction force vector using Ht map: f3 = Ht^T * fric
285  f3(1:edof) = matmul( transpose(ht(1:2,1:edof)), fric(1:2) )
286 
287  ! Distribute friction force (negative for equilibrium)
288  cttforce(1:edof) = -f3(1:edof)
289 
290  end subroutine computefrictionforce_alag
291 
293  subroutine gettrialfricforceandcheckfricstate(ctstate,fcoeff,tPenalty,lagrange)
294 
295  type(tcontactstate) :: ctstate
296  real(kind=kreal) :: fcoeff, tpenalty
297  real(kind=kreal) :: lagrange
298 
299  integer(kind=kint) :: i
300  real(kind=kreal) :: tf_yield
301 
302  do i = 1, 3
303  ctstate%tangentForce_trial(i) = ctstate%tangentForce1(i) + tpenalty*ctstate%reldisp(i)
304  enddo
305 
306  tf_yield = fcoeff*dabs(lagrange)
307  if(ctstate%state == contactslip) tf_yield =0.99d0*tf_yield
308  if( dsqrt(dot_product(ctstate%tangentForce_trial,ctstate%tangentForce_trial)) <= tf_yield ) then
309  ctstate%state = contactstick
310  else
311  ctstate%state = contactslip
312  endif
313 
315 
317  subroutine update_tangentforce(etype,nn,elemt0,elemt,cstate)
318  integer, intent(in) :: etype
319  integer, intent(in) :: nn
320  real(kind=kreal),intent(in) :: elemt0(3,nn)
321  real(kind=kreal),intent(in) :: elemt(3,nn)
322  type(tcontactstate), intent(inout) :: cstate
323 
324  integer :: i
325  real(kind=kreal) :: tangent0(3,2), tangent(3,2) ! base vectors in tangent space
326  real(kind=kreal) :: coeff(2), norm, norm_tmp
327  real(kind=kreal) :: tangentforce_tmp(3)
328 
329  call tangentbase( etype, nn, cstate%lpos(1:2), elemt0, tangent0 )
330  call tangentbase( etype, nn, cstate%lpos(1:2), elemt, tangent )
331 
332  !project tangentforce to base vector tangent0
333  do i=1,2
334  coeff(i) = dot_product(cstate%tangentForce(1:3),tangent0(1:3,i))
335  coeff(i) = coeff(i)/dot_product(tangent0(1:3,i),tangent0(1:3,i))
336  enddo
337  tangentforce_tmp(1:3) = coeff(1)*tangent0(1:3,1) + coeff(2)*tangent0(1:3,2)
338  norm_tmp = dsqrt(dot_product(tangentforce_tmp,tangentforce_tmp))
339  !adjust tangent force of slave point which moved over element boundary
340  if( norm_tmp > 1.d-6 ) then
341  norm = dsqrt(dot_product(cstate%tangentForce,cstate%tangentForce))
342  coeff(1:2) = (norm/norm_tmp)*coeff(1:2)
343  end if
344 
345  !set rotated tangentforce to tangentforce1
346  cstate%tangentForce1(1:3) = coeff(1)*tangent(1:3,1) + coeff(2)*tangent(1:3,2)
347 
348  end subroutine update_tangentforce
349 
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
subroutine getshapefunc(fetype, localcoord, func)
Calculate the shape function in natural coordinate system.
Definition: element.f90:647
subroutine tangentbase(fetype, nn, localcoord, elecoord, tangent)
Calculate base vector of tangent space of 3d surface.
Definition: element.f90:933
Definition: hecmw.f90:6
Common utilities for contact element calculations.
subroutine, public computetm_tt(ctState, tSurf, fcoeff, Tm, Tt)
Compute Tm (relative displacement mapping) and optionally Tt (tangential mapping) This subroutine con...
subroutine, public computecontactmaps_alag(ctState, tSurf, ele, Bn, metric, Ht, Gt)
Compute contact maps for ALag method (normal distribution and tangential displacement maps) This subr...
subroutine, public gettrialfricforceandcheckfricstate(ctstate, fcoeff, tPenalty, lagrange)
This subroutine calculates trial friction force and checks friction state.
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,...
subroutine, public update_tangentforce(etype, nn, elemt0, elemt, cstate)
This subroutine find the projection of a slave point onto master surface.
This module provides geometric calculations for contact.
subroutine dispincrematrix(pos, etype, nnode, ele, tangent, tensor, matrix)
This subroutine calculate the relation between global disp and displacement along natural coordinate ...
This module manages the data structure for contact calculation.
integer, parameter contactslip
integer, parameter contactstick
This structure records contact status.