FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_contact_elem_alag.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
10  use msurfelement
13  implicit none
14 
15  public :: getcontactstiffness_alag
17  public :: gettiedstiffness_alag
18  public :: gettiednodalforce_alag
20 
21 contains
22 
23  subroutine getcontactstiffness_alag(cstate, tSurf, ele, mu, mut, fcoeff, symm, stiff, force)
24 
25  type(tcontactstate), intent(inout) :: cstate
26  type(tsurfelement), intent(in) :: tsurf
27  real(kind=kreal), intent(in) :: ele(:,:)
28  real(kind=kreal), intent(in) :: mu, mut
29  real(kind=kreal), intent(in) :: fcoeff
30  logical, intent(in) :: symm
31  real(kind=kreal), intent(out) :: stiff(:,:)
32  real(kind=kreal), intent(out) :: force(:)
33 
34  integer :: i, j, nnode
35  real(kind=kreal) :: bn(size(tsurf%nodes)*3+3), ht(2,size(tsurf%nodes)*3+3), gt(2,size(tsurf%nodes)*3+3)
36  real(kind=kreal) :: metric(2,2)
37  real(kind=kreal) :: a(2,2)
38  real(kind=kreal) :: norm_trial, alpha_proj, that(2)
39  real(kind=kreal) :: k_fric(size(tsurf%nodes)*3+3,size(tsurf%nodes)*3+3)
40  real(kind=kreal) :: hta(2,size(tsurf%nodes)*3+3), tmp_vec(2)
41  real(kind=kreal) :: dummy_force(size(tsurf%nodes)*3+3)
42  real(kind=kreal) :: zero_disp(size(tsurf%nodes)*3+3)
43 
44  nnode = size(tsurf%nodes)
45 
46  ! Use common mapping routine to compute Bn, metric, Ht, Gt
47  call computecontactmaps_alag(cstate, tsurf, ele, bn, metric, ht, gt)
48 
49  ! Normal stiffness: stiff = mu * Bn * Bn^T
50  do j = 1, nnode*3+3
51  do i = 1, nnode*3+3
52  stiff(i,j) = mu * bn(i) * bn(j)
53  enddo
54  enddo
55  force(1:nnode*3+3) = bn(:)
56 
57  ! frictional component
58  if( fcoeff /= 0.d0 ) then
59  ! Compute trial friction info for consistent tangent
60  ! Use zero displacement to get projection parameters from current multiplier state
61  zero_disp = 0.0d0 ! zero displacement increment
62  call computefrictionforce_alag(cstate, fcoeff, cstate%multiplier(1), metric, &
63  ht, gt, zero_disp, nnode*3+3, dummy_force, &
64  mut, norm_trial=norm_trial, alpha=alpha_proj, that=that)
65 
66  ! Construct 2D local tangent operator A
67  if( cstate%multiplier(1) <= 0.0d0 .or. alpha_proj <= 1.0d-20 ) then
68  ! No friction: A = 0
69  a = 0.0d0
70  else if( alpha_proj >= 0.999d0 .or. norm_trial < 1.0d-20 ) then
71  ! Stick: A = mut*I
72  a(1,1) = mut
73  a(2,2) = mut
74  a(1,2) = 0.0d0
75  a(2,1) = 0.0d0
76  else
77  ! Slip: A = alpha*mut*(I - that⊗that)
78  a(1,1) = alpha_proj * mut * (1.0d0 - that(1)*that(1))
79  a(1,2) = alpha_proj * mut * (-that(1)*that(2))
80  a(2,1) = alpha_proj * mut * (-that(2)*that(1))
81  a(2,2) = alpha_proj * mut * (1.0d0 - that(2)*that(2))
82  endif
83 
84  ! Compute friction stiffness: K_fric = Ht^T * A * Gt
85  ! First: HtA = A^T * Ht (2 x edof)
86  hta(1:2,1:nnode*3+3) = matmul(transpose(a), ht(1:2,1:nnode*3+3))
87  ! Then: K_fric = Ht^T * HtA = Ht^T * A^T * Ht (but we want Ht^T * A * Gt)
88  ! Correct: K_fric(i,j) = Ht(:,i)^T * A * Gt(:,j)
89  do j = 1, nnode*3+3
90  tmp_vec = matmul(a, gt(1:2,j))
91  do i = 1, nnode*3+3
92  k_fric(i,j) = dot_product(ht(1:2,i), tmp_vec)
93  enddo
94  enddo
95 
96  stiff(1:nnode*3+3,1:nnode*3+3) = stiff(1:nnode*3+3,1:nnode*3+3) + k_fric(1:nnode*3+3,1:nnode*3+3)
97  endif
98 
99  end subroutine getcontactstiffness_alag
100 
101  subroutine getcontactnodalforce_alag(ctState,tSurf,ndCoord,ndDu,mu,mut,fcoeff,lagrange,ctNForce,ctTForce,cflag)
102 
103  use msurfelement
104  type(tcontactstate) :: ctstate
105  type(tsurfelement) :: tsurf
106  integer(kind=kint) :: nnode
107  integer(kind=kint) :: j
108  real(kind=kreal), intent(in) :: mu, mut
109  real(kind=kreal) :: fcoeff
110  real(kind=kreal) :: lagrange
111  real(kind=kreal) :: ndcoord(:), nddu(:)
112  real(kind=kreal) :: ctnforce(:)
113  real(kind=kreal) :: cttforce(:)
114  logical :: cflag
115 
116  real(kind=kreal) :: normal(3)
117  real(kind=kreal) :: bn(3*l_max_elem_node+3)
118  real(kind=kreal) :: ht(2,3*l_max_elem_node+3), gt(2,3*l_max_elem_node+3)
119  real(kind=kreal) :: elemcrd(3, l_max_elem_node)
120  real(kind=kreal) :: edisp(3*l_max_elem_node+3)
121  real(kind=kreal) :: dgn, nrlforce
122  real(kind=kreal) :: metric(2,2)
123  integer(kind=kint) :: edof
124 
125  nnode = size(tsurf%nodes)
126  edof = nnode*3+3
127 
128  ctnforce = 0.0d0
129  cttforce = 0.0d0
130 
131  normal(1:3) = ctstate%direction(1:3)
132 
133  ! Prepare elemcrd = ndCoord - ndDu (i.e., coord + disp) for computeContactMaps_ALag
134  do j = 1, nnode
135  elemcrd(1:3, j) = ndcoord(j*3+1:j*3+3) - nddu(j*3+1:j*3+3)
136  enddo
137 
138  ! Use common mapping routine to compute Bn, metric, Ht, Gt
139  call computecontactmaps_alag(ctstate, tsurf, elemcrd(:,1:nnode), &
140  bn, metric, ht, gt)
141 
142  ! Normal gap: dgn = Bn^T * ndCoord (using normal distribution vector)
143  dgn = dot_product( bn(1:edof), ndcoord(1:edof) )
144 
145  ! Normal force: multiplier + penalty * gap
146  nrlforce = ctstate%multiplier(1) + mu*dgn
147 
148  ! Distribute normal force using Bn: ctNForce = -nrlforce * Bn
149  ctnforce(1:edof) = -nrlforce * bn(1:edof)
150 
151  ! Lagrange row (not used in ALagrange, set to 0)
152  ctnforce((nnode+1)*3+1) = 0.d0
153 
154  if( fcoeff == 0.d0 ) return
155 
156  ! --- Tangent component ---
157 
158  ! Prepare edisp from ndDu
159  edisp(1:3) = nddu(1:3) ! slave
160  do j = 1, nnode
161  edisp(j*3+1:j*3+3) = nddu(j*3+1:j*3+3) ! master nodes
162  enddo
163 
164  ! Compute friction force using common routine
165  call computefrictionforce_alag(ctstate, fcoeff, ctstate%multiplier(1), metric, &
166  ht, gt, edisp, edof, cttforce, &
167  mut)
168 
169  ! Lagrange row (not used in ALagrange, set to 0)
170  cttforce((nnode+1)*3+1) = 0.d0
171 
172  end subroutine getcontactnodalforce_alag
173 
174  subroutine updatecontactmultiplier_alag(ctState,ndLocal,coord,disp,ddisp,&
175  & mu,mut,fcoeff,etype,lgnt,ctchanged,ctNForce,ctTForce,jump_ratio)
176 
177  type(tcontactstate), intent(inout) :: ctstate
178  integer(kind=kint), intent(in) :: ndlocal(:)
179  real(kind=kreal), intent(in) :: coord(:)
180  real(kind=kreal), intent(in) :: disp(:)
181  real(kind=kreal), intent(in) :: ddisp(:)
182  real(kind=kreal), intent(in) :: mu, mut
183  real(kind=kreal), intent(in) :: fcoeff
184  integer(kind=kint), intent(in) :: etype
185  real(kind=kreal), intent(inout) :: lgnt(2)
186  logical, intent(inout) :: ctchanged
187  real(kind=kreal), intent(out) :: ctnforce(:)
188  real(kind=kreal), intent(out) :: cttforce(:)
189  real(kind=kreal), intent(out) :: jump_ratio
190 
191  integer(kind=kint) :: nnode
192  integer(kind=kint) :: slave, j
193  real(kind=kreal) :: bn(3*l_max_elem_node+3)
194  real(kind=kreal) :: ht(2,3*l_max_elem_node+3), gt(2,3*l_max_elem_node+3)
195  real(kind=kreal) :: elemcrd(3,l_max_elem_node)
196  real(kind=kreal) :: curpos(3*l_max_elem_node+3)
197  real(kind=kreal) :: edisp(3*l_max_elem_node+3)
198  real(kind=kreal) :: dgn, nrlforce
199  real(kind=kreal) :: metric(2,2)
200  real(kind=kreal) :: dxy(2)
201  integer(kind=kint) :: edof
202  type(tsurfelement) :: tsurf_tmp
203 
204  nnode = size(ndlocal) - 1
205  slave = ndlocal(1)
206  edof = nnode*3+3
207 
208  ctnforce = 0.0d0
209  cttforce = 0.0d0
210 
211  ! Prepare elemcrd (coord+disp) and current positions (coord+disp+ddisp)
212  curpos(1:3) = coord(3*slave-2:3*slave) + disp(3*slave-2:3*slave) + ddisp(3*slave-2:3*slave)
213  edisp(1:3) = ddisp(3*slave-2:3*slave)
214  do j = 1, nnode
215  elemcrd(1:3,j) = coord(3*ndlocal(j+1)-2:3*ndlocal(j+1)) + disp(3*ndlocal(j+1)-2:3*ndlocal(j+1))
216  curpos(j*3+1:j*3+3) = elemcrd(1:3,j) + ddisp(3*ndlocal(j+1)-2:3*ndlocal(j+1))
217  edisp(j*3+1:j*3+3) = ddisp(3*ndlocal(j+1)-2:3*ndlocal(j+1))
218  enddo
219 
220  ! Build temporary tSurf for computeContactMaps_ALag
221  allocate(tsurf_tmp%nodes(nnode))
222  tsurf_tmp%etype = etype
223  tsurf_tmp%nodes = ndlocal(2:nnode+1)
224 
225  ! Use common mapping routine to compute Bn, metric, Ht, Gt
226  call computecontactmaps_alag(ctstate, tsurf_tmp, elemcrd(:,1:nnode), &
227  bn, metric, ht, gt)
228 
229  ! Normal gap: dgn = Bn^T * curpos (using normal distribution vector)
230  dgn = dot_product( bn(1:edof), curpos(1:edof) )
231 
232  ! Update multiplier and working distance
233  ctstate%wkdist = -dgn
234  ctstate%multiplier(1) = ctstate%multiplier(1) - mu*ctstate%wkdist
235  ctstate%distance = ctstate%wkdist
236  lgnt(1) = lgnt(1) - ctstate%wkdist
237 
238  ! Normal force: use updated multiplier
239  nrlforce = ctstate%multiplier(1)
240 
241  ! Distribute normal force using Bn: ctNForce = -nrlforce * Bn
242  ctnforce(1:edof) = -nrlforce * bn(1:edof)
243 
244  if( fcoeff == 0.d0 ) return
245 
246  ! --- Tangent component ---
247 
248  ! Compute friction force with multiplier update and state check
249  call computefrictionforce_alag(ctstate, fcoeff, ctstate%multiplier(1), metric, &
250  ht, gt, edisp, edof, cttforce, &
251  mut, &
252  update_multiplier=.true., slave_id=slave, ctchanged=ctchanged, &
253  jump_ratio=jump_ratio)
254 
255  ! Tangent displacement for convergence check: use Gt to project curpos directly
256  dxy = matmul( gt(:,1:edof), curpos(1:edof) )
257  lgnt(2) = lgnt(2) + dsqrt( dxy(1)*dxy(1) + dxy(2)*dxy(2) )
258 
259  ! Clean up temporary surface element
260  deallocate(tsurf_tmp%nodes)
261 
262  end subroutine updatecontactmultiplier_alag
263 
264  subroutine gettiedstiffness_alag(cstate, tSurf, mu, stiff, force)
265 
266  type(tcontactstate), intent(in) :: cstate
267  type(tsurfelement), intent(in) :: tsurf
268  real(kind=kreal), intent(in) :: mu
269  real(kind=kreal), intent(out) :: stiff(:,:)
270  real(kind=kreal), intent(out) :: force(:)
271 
272  integer :: i, j, nnode, edof
273  real(kind=kreal) :: tm(3, 3*(l_max_surface_node+1))
274  real(kind=kreal) :: tt(3, 3*(l_max_surface_node+1))
275 
276  nnode = size(tsurf%nodes)
277  edof = nnode*3+3
278 
279  stiff = 0.d0
280 
281  ! Use common mapping routine to compute Tm
282  call computetm_tt(cstate, tsurf, 0.0d0, tm, tt)
283 
284  ! Tied stiffness: stiff = mu * Tm^T * Tm
285  do j = 1, edof
286  do i = 1, edof
287  stiff(i,j) = mu * dot_product(tm(1:3,i), tm(1:3,j))
288  enddo
289  enddo
290  force(1:edof) = 0.d0 ! not used for tied (3-direction constraint)
291 
292  end subroutine gettiedstiffness_alag
293 
294  subroutine gettiednodalforce_alag(ctState,tSurf,ndu,mu,ctNForce,ctTForce)
295 
296  use msurfelement
297  type(tcontactstate) :: ctstate
298  type(tsurfelement) :: tsurf
299  integer(kind=kint) :: nnode
300  real(kind=kreal) :: ndu(:)
301  real(kind=kreal), intent(in) :: mu
302  real(kind=kreal) :: ctnforce(:)
303  real(kind=kreal) :: cttforce(:)
304 
305  integer(kind=kint) :: edof
306  real(kind=kreal) :: tm(3, 3*(l_max_surface_node+1))
307  real(kind=kreal) :: tt(3, 3*(l_max_surface_node+1))
308  real(kind=kreal) :: dg(3)
309  real(kind=kreal) :: nrlforce(3)
310 
311  nnode = size(tsurf%nodes)
312  edof = nnode*3+3
313 
314  ctnforce = 0.0d0
315  cttforce = 0.0d0
316 
317  ! Use common mapping routine to compute Tm
318  call computetm_tt(ctstate, tsurf, 0.0d0, tm, tt)
319 
320  ! Gap vector: dg = Tm * ndu (3-component relative displacement)
321  dg(1:3) = matmul(tm(1:3, 1:edof), ndu(1:edof))
322 
323  ! Force: multiplier + penalty * gap (3 components)
324  nrlforce(1:3) = ctstate%multiplier(1:3) + mu*dg(1:3)
325 
326  ! Distribute force: ctNForce = -Tm^T * nrlforce
327  ctnforce(1:edof) = -matmul(transpose(tm(1:3, 1:edof)), nrlforce(1:3))
328 
329  ! Lagrange row (not used in ALagrange, set to 0)
330  ctnforce((nnode+1)*3+1) = 0.d0
331  cttforce((nnode+1)*3+1) = 0.d0
332 
333  end subroutine gettiednodalforce_alag
334 
335 end module m_fstr_contact_elem_alag
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
Definition: hecmw.f90:6
Alag method implementations for contact element calculations.
subroutine, public updatecontactmultiplier_alag(ctState, ndLocal, coord, disp, ddisp, mu, mut, fcoeff, etype, lgnt, ctchanged, ctNForce, ctTForce, jump_ratio)
subroutine, public gettiedstiffness_alag(cstate, tSurf, mu, stiff, force)
subroutine, public getcontactnodalforce_alag(ctState, tSurf, ndCoord, ndDu, mu, mut, fcoeff, lagrange, ctNForce, ctTForce, cflag)
subroutine, public gettiednodalforce_alag(ctState, tSurf, ndu, mu, ctNForce, ctTForce)
subroutine, public getcontactstiffness_alag(cstate, tSurf, ele, mu, mut, fcoeff, symm, stiff, force)
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 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,...
This module provides geometric calculations for contact.
This module manages the data structure for contact calculation.
This module manages surface elements in 3D It provides basic definition of surface elements (triangla...
Definition: surf_ele.f90:8
integer(kind=kint), parameter l_max_elem_node
Definition: surf_ele.f90:17
integer(kind=kint), parameter l_max_surface_node
Definition: surf_ele.f90:16
This structure records contact status.
Structure to define surface group.
Definition: surf_ele.f90:23