FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_contact_elem_slag.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
14  implicit none
15 
16  public :: getcontactstiffness_slag
18  public :: gettiedstiffness_slag
19  public :: gettiednodalforce_slag
20 
21 contains
22 
23  subroutine getcontactstiffness_slag(ctState,tSurf,iter,tPenalty,fcoeff,lagrange,stiffness,smoothing_type)
24 
25  use msurfelement
26  type(tcontactstate) :: ctstate
27  type(tsurfelement) :: tsurf
28  integer(kind=kint) :: iter
29  integer(kind=kint) :: nnode
30  integer(kind=kint) :: i, j
31  real(kind=kreal) :: normal(3)
32  real(kind=kreal) :: ntm((l_max_surface_node + 1)*3)
33  real(kind=kreal) :: fcoeff, tpenalty
34  real(kind=kreal) :: lagrange
35  real(kind=kreal) :: tf_trial(3), length_tft
36  real(kind=kreal) :: tangent(3)
37  real(kind=kreal) :: stiffness(:,:)
38  real(kind=kreal) :: tm(3, 3*(l_max_surface_node+1)), tt(3, 3*(l_max_surface_node+1))
39  real(kind=kreal) :: q((l_max_surface_node+1)*3)
40  real(kind=kreal) :: dot_tn
41  integer(kind=kint), optional :: smoothing_type
42 
43  nnode = size(tsurf%nodes)
44 
45  stiffness = 0.0d0
46 
47  ! Compute Tm and Tt matrices using standard shape functions
48  call computetm_tt(ctstate, tsurf, fcoeff, tm, tt, smoothing_type)
49 
50  normal(1:3) = ctstate%direction(1:3)
51 
52  ! Compute nTm = Tm^T * n for KKT constraint (Lagrange multiplier row/column)
53  ntm(1:(nnode+1)*3) = matmul(transpose(tm(1:3, 1:3*(nnode+1))), normal)
54 
55  i = (nnode+1)*3 + 1
56  do j = 1, (nnode+1)*3
57  stiffness(i,j) = ntm(j); stiffness(j,i) = ntm(j)
58  enddo
59 
60  if( fcoeff /= 0.0d0 ) then
61  if( lagrange>0.0d0 .or. iter==1 ) then
62 
63  ! Stick friction: K_uu^stick = ρ * Tt^T * Tt
64  do i = 1, (nnode+1)*3
65  do j = 1, i
66  stiffness(i,j) = stiffness(i,j) + tpenalty * dot_product(tt(1:3,i), tt(1:3,j))
67  ! Fill symmetric part
68  if (i /= j) then
69  stiffness(j,i) = stiffness(i,j)
70  endif
71  enddo
72  enddo
73 
74  if( ctstate%state == contactslip ) then
75 
76  tf_trial(1:3) = ctstate%tangentForce_trial(1:3)
77  length_tft = dsqrt(dot_product(tf_trial,tf_trial))
78  tangent(1:3) = tf_trial(1:3)/length_tft
79 
80  ! Compute q = Tm^T * tangent
81  q(1:(nnode+1)*3) = matmul(transpose(tm(1:3, 1:3*(nnode+1))), tangent)
82 
83  ! Slip correction:
84  dot_tn = dot_product(tangent, normal)
85 
86  do i = 1, (nnode+1)*3
87  do j = 1, i
88  stiffness(i,j) = stiffness(i,j) + tpenalty * (-q(i)*q(j) + q(i)*ntm(j)*dot_tn)
89  if (i /= j) then
90  stiffness(j,i) = stiffness(i,j)
91  endif
92  enddo
93  enddo
94 
95  ! Scale (u,u) block by (μλ/|t_tr|)
96  stiffness(1:(nnode+1)*3,1:(nnode+1)*3) = (fcoeff*lagrange/length_tft) * stiffness(1:(nnode+1)*3,1:(nnode+1)*3)
97 
98  ! (u,λ) cross term:)
99  stiffness(1:(nnode+1)*3, (nnode+1)*3+1) = stiffness(1:(nnode+1)*3, (nnode+1)*3+1) + fcoeff * q(1:(nnode+1)*3)
100 
101  endif
102  endif
103  endif
104 
105  end subroutine getcontactstiffness_slag
106 
107  subroutine getcontactnodalforce_slag(ctState,tSurf,ndCoord,ndDu,tPenalty,fcoeff,lagrange,ctNForce,ctTForce,cflag,smoothing_type)
108 
109  use msurfelement
110  type(tcontactstate) :: ctstate
111  type(tsurfelement) :: tsurf
112  integer(kind=kint) :: nnode
113  integer(kind=kint), optional :: smoothing_type
114  integer(kind=kint) :: j
115  real(kind=kreal) :: fcoeff, tpenalty
116  real(kind=kreal) :: lagrange
117  real(kind=kreal) :: ndcoord(:), nddu(:)
118  real(kind=kreal) :: normal(3)
119  real(kind=kreal) :: ntm((l_max_surface_node + 1)*3)
120  real(kind=kreal) :: tf_trial(3), length_tft, tangent(3), tf_final(3)
121  real(kind=kreal) :: ctnforce(:)
122  real(kind=kreal) :: cttforce(:)
123  logical :: cflag
124 
125  real(kind=kreal) :: tm(3, 3*(l_max_surface_node+1)), tt(3, 3*(l_max_surface_node+1))
126 
127  ctstate%multiplier(1) = lagrange
128 
129  nnode = size(tsurf%nodes)
130 
131  ctnforce = 0.0d0
132  cttforce = 0.0d0
133 
134  ! Compute Tm matrix
135  call computetm_tt(ctstate, tsurf, fcoeff, tm, tt, smoothing_type)
136 
137  normal(1:3) = ctstate%direction(1:3)
138 
139  ! Compute nTm = Tm^T * normal
140  ntm(1:(nnode+1)*3) = -matmul(transpose(tm(1:3, 1:3*(nnode+1))), normal)
141 
142  do j = 1, (nnode+1)*3
143  ctnforce(j) = lagrange*ntm(j)
144  enddo
145  j = (nnode+1)*3 + 1
146  ctnforce(j) = dot_product(ntm(1:(nnode+1)*3),ndcoord(1:(nnode+1)*3))
147 
148  if(fcoeff /= 0.0d0 .and. lagrange > 0.0d0)then
149 
150  if( cflag ) then !calc tf_final and set it to tangentForce_final
151  ctstate%reldisp(1:3) = matmul(tt(1:3,1:3*(nnode+1)), nddu(1:3*(nnode+1)))
152 
153  call gettrialfricforceandcheckfricstate(ctstate,fcoeff,tpenalty,lagrange)
154 
155  if( ctstate%state == contactstick ) then
156  tf_final(1:3) = ctstate%tangentForce_trial(1:3)
157  elseif( ctstate%state == contactslip ) then
158  tf_trial(1:3) = ctstate%tangentForce_trial(1:3)
159  length_tft = dsqrt(dot_product(tf_trial,tf_trial))
160  tangent(1:3) = tf_trial(1:3)/length_tft
161  tf_final(1:3) = fcoeff*dabs(lagrange)*tangent(1:3)
162  endif
163  ctstate%tangentForce_final(1:3) = tf_final(1:3)
164  else ! just set tangentForce_final to tf_final (used for fstr_ass_load_contact)
165  tf_final(1:3) = ctstate%tangentForce_final(1:3)
166  endif
167 
168  ! Distribute friction force using Tm: ctTForce = -Tm^T * tf_final
169  cttforce(1:(nnode+1)*3) = -matmul(transpose(tm(1:3, 1:3*(nnode+1))), tf_final)
170 
171  endif
172 
173  end subroutine getcontactnodalforce_slag
174 
175  subroutine gettiedstiffness_slag(ctState,tSurf,idof,stiffness,smoothing_type)
176  use msurfelement
177  type(tcontactstate) :: ctstate
178  type(tsurfelement) :: tsurf
179  integer(kind=kint) :: nnode
180  integer(kind=kint) :: idof
181  real(kind=kreal) :: stiffness(:,:)
182  integer(kind=kint), optional :: smoothing_type
183 
184  integer(kind=kint) :: i, j
185  real(kind=kreal) :: ntm((l_max_surface_node+1)*3)
186  real(kind=kreal) :: tm(3, 3*(l_max_surface_node+1)), tt(3, 3*(l_max_surface_node+1))
187  real(kind=kreal) :: e_idof(3)
188 
189  nnode = size(tsurf%nodes)
190  stiffness = 0.0d0
191 
192  ! Compute Tm matrix
193  call computetm_tt(ctstate, tsurf, 0.0d0, tm, tt, smoothing_type)
194 
195  ! Create unit vector in idof direction
196  e_idof = 0.0d0
197  e_idof(idof) = 1.0d0
198 
199  ! Compute nTm = Tm^T * e_idof
200  ! This gives constraint direction for tied contact
201  ntm(1:(nnode+1)*3) = matmul(transpose(tm(1:3, 1:3*(nnode+1))), e_idof)
202 
203  i = (nnode+1)*3 + 1
204  do j = 1, (nnode+1)*3
205  stiffness(i,j) = ntm(j); stiffness(j,i) = ntm(j)
206  enddo
207 
208  end subroutine gettiedstiffness_slag
209 
210  subroutine gettiednodalforce_slag(ctState,tSurf,idof,ndu,lagrange,ctNForce,ctTForce,smoothing_type)
211  use msurfelement
212  type(tcontactstate) :: ctstate
213  type(tsurfelement) :: tsurf
214  integer(kind=kint) :: nnode
215  integer(kind=kint) :: idof
216  real(kind=kreal) :: ndu(:)
217  real(kind=kreal) :: lagrange
218  real(kind=kreal) :: ctnforce(:)
219  real(kind=kreal) :: cttforce(:)
220  integer(kind=kint), optional :: smoothing_type
221 
222  integer(kind=kint) :: j
223  real(kind=kreal) :: normal(3)
224  real(kind=kreal) :: ntm((l_max_surface_node + 1)*3)
225  real(kind=kreal) :: ttm((l_max_surface_node + 1)*3)
226  real(kind=kreal) :: tm(3, 3*(l_max_surface_node+1)), tt(3, 3*(l_max_surface_node+1))
227  real(kind=kreal) :: e_idof(3), e_normal(3), e_tangent(3)
228 
229  nnode = size(tsurf%nodes)
230 
231  ctnforce = 0.0d0
232  cttforce = 0.0d0
233 
234  ! Compute Tm matrix
235  call computetm_tt(ctstate, tsurf, 0.0d0, tm, tt, smoothing_type)
236 
237  normal(1:3) = ctstate%direction(1:3)
238 
239  ! Create unit vector in idof direction
240  e_idof = 0.0d0
241  e_idof(idof) = 1.0d0
242 
243  ! Normal component: -normal(idof) * normal (OLD: nTm(1:3) = -normal(idof)*normal(1:3))
244  e_normal = -normal(idof) * normal
245 
246  ! Tangential component: -e_idof - e_normal (OLD: tTm(idof)=-1.0; tTm(1:3)-=nTm(1:3))
247  e_tangent = -e_idof - e_normal
248 
249  ! Compute nTm and tTm using Tm^T
250  ntm(1:(nnode+1)*3) = matmul(transpose(tm(1:3, 1:3*(nnode+1))), e_normal)
251  ttm(1:(nnode+1)*3) = matmul(transpose(tm(1:3, 1:3*(nnode+1))), e_tangent)
252 
253  do j = 1, (nnode+1)*3
254  ctnforce(j) = lagrange*ntm(j)
255  cttforce(j) = lagrange*ttm(j)
256  enddo
257  j = (nnode+1)*3 + 1
258  ctnforce(j) = dot_product(ntm(1:(nnode+1)*3),ndu(1:(nnode+1)*3))
259  cttforce(j) = dot_product(ttm(1:(nnode+1)*3),ndu(1:(nnode+1)*3))
260  end subroutine gettiednodalforce_slag
261 
262 end module m_fstr_contact_elem_slag
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
Definition: hecmw.f90:6
Common utilities for contact element calculations.
subroutine, public gettrialfricforceandcheckfricstate(ctstate, fcoeff, tPenalty, lagrange)
This subroutine calculates trial friction force and checks friction state.
subroutine, public computetm_tt(ctState, tSurf, fcoeff, Tm, Tt, smoothing_type, Bn)
Compute Tm (relative displacement mapping) and optionally Tt (tangential mapping) This subroutine con...
Slag method implementations for contact element calculations.
subroutine, public gettiednodalforce_slag(ctState, tSurf, idof, ndu, lagrange, ctNForce, ctTForce, smoothing_type)
subroutine, public getcontactnodalforce_slag(ctState, tSurf, ndCoord, ndDu, tPenalty, fcoeff, lagrange, ctNForce, ctTForce, cflag, smoothing_type)
subroutine, public gettiedstiffness_slag(ctState, tSurf, idof, stiffness, smoothing_type)
subroutine, public getcontactstiffness_slag(ctState, tSurf, iter, tPenalty, fcoeff, lagrange, stiffness, smoothing_type)
This module provides geometric calculations for contact.
Contact surface smoothing using Nagata patch interpolation.
This module manages the data structure for contact calculation.
integer, parameter contactslip
integer, parameter contactstick
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_surface_node
Definition: surf_ele.f90:16
This structure records contact status.
Structure to define surface group.
Definition: surf_ele.f90:23