FrontISTR  5.8.0
Large-scale structural analysis program with finit element method
fstr_contact_def.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 !-------------------------------------------------------------------------------
13 
14  use hecmw
15  use msurfelement
16  use m_contact_lib
18  use bucket_search
19  use mcontactparam
21 
22  implicit none
23 
25  type tcontact
26  ! following contact definition
27  character(len=HECMW_NAME_LEN) :: name
28  integer :: ctype
29  integer :: group
30  character(len=HECMW_NAME_LEN) :: pair_name
31  integer :: surf_id1, surf_id2
32  integer :: surf_id1_sgrp
33  type(tsurfelement), pointer :: master(:)=>null()
34  integer, pointer :: slave(:)=>null()
35  real(kind=kreal) :: fcoeff
36  real(kind=kreal) :: tpenalty
37 
38  real(kind=kreal) :: ctime
39  integer(kind=kint) :: if_type
40  real(kind=kreal) :: if_etime
41  real(kind=kreal) :: initial_pos
42  real(kind=kreal) :: end_pos
43  ! following algorithm
44  ! -1: not initialized
45  ! 1: TIED-Just rigidly fixed the two surfaces
46  ! 2: GLUED-Distance between the two surfaces to zero and glue them
47  ! 3: SSLID-Small sliding contact( no position but with contact state change)
48  ! 4: FSLID-Finite sliding contact (both changes in contact state and position possible)
49  integer :: algtype
50 
51  logical :: mpced
52  logical :: symmetric
53 
54  ! following contact state
55  type(tcontactstate), pointer :: states(:)=>null()
56 
57  type(hecmwst_contact_comm) :: comm
58  type(bucketdb) :: master_bktdb
59 
60  type(tcontactparam), pointer :: cparam=>null()
61  end type tcontact
62 
64  logical :: active
65  integer(kind=kint) :: contact2free
66  integer(kind=kint) :: contact2neighbor
67  integer(kind=kint) :: contact2difflpos
68  integer(kind=kint) :: free2contact
69  integer(kind=kint) :: contactnode_previous
70  integer(kind=kint) :: contactnode_current
72 
73  private :: is_mpc_available
74  private :: is_active_contact
75  private :: cal_node_normal
76  private :: track_contact_position
77  private :: reset_contact_force
78 
79 contains
80 
82  logical function is_mpc_available( contact )
83  type(tcontact), intent(in) :: contact
84  is_mpc_available = .true.
85  if( contact%fcoeff/=0.d0 ) is_mpc_available = .false.
86  end function
87 
89  subroutine fstr_write_contact( file, contact )
90  integer(kind=kint), intent(in) :: file
91  type(tcontact), intent(in) :: contact
92  integer :: i
93  write(file,*) "CONTACT:", contact%ctype,contact%group,trim(contact%pair_name),contact%fcoeff
94  write(file,*) "---Slave----"
95  write(file,*) "num.slave",size(contact%slave)
96  if( associated(contact%slave) ) then
97  do i=1,size(contact%slave)
98  write(file, *) contact%slave(i)
99  enddo
100  endif
101  write(file,*) "----master---"
102  write(file,*) "num.master",size(contact%master)
103  if( associated(contact%master) ) then
104  do i=1,size(contact%master)
105  call write_surf( file, contact%master(i) )
106  enddo
107  endif
108  end subroutine
109 
111  subroutine fstr_contact_finalize( contact )
112  type(tcontact), intent(inout) :: contact
113  integer :: i
114  if( associated( contact%slave ) ) deallocate(contact%slave)
115  if( associated( contact%master ) ) then
116  do i=1,size( contact%master )
117  call finalize_surf( contact%master(i) )
118  enddo
119  deallocate(contact%master)
120  endif
121  if( associated(contact%states) ) deallocate(contact%states)
122  call hecmw_contact_comm_finalize(contact%comm)
123  call bucketdb_finalize( contact%master_bktDB )
124  end subroutine
125 
127  logical function fstr_contact_check( contact, hecMESH )
128  type(tcontact), intent(inout) :: contact
129  type(hecmwst_local_mesh), pointer :: hecmesh
130 
131  integer :: i
132  logical :: isfind
133 
134  fstr_contact_check = .false.
135 
136  ! if contact pair exist?
137  isfind = .false.
138  do i=1,hecmesh%contact_pair%n_pair
139  if( hecmesh%contact_pair%name(i) == contact%pair_name ) then
140  contact%ctype = hecmesh%contact_pair%type(i)
141  contact%surf_id1 = hecmesh%contact_pair%slave_grp_id(i)
142  contact%surf_id2 = hecmesh%contact_pair%master_grp_id(i)
143  contact%surf_id1_sgrp = hecmesh%contact_pair%slave_orisgrp_id(i)
144  isfind = .true.
145  endif
146  enddo
147  if( .not. isfind ) return;
148  if( contact%fcoeff<=0.d0 ) contact%fcoeff=0.d0
149  if( contact%ctype < 1 .and. contact%ctype > 3 ) return
150  if( contact%group<=0 ) return
151 
152  fstr_contact_check = .true.
153  end function
154 
156  logical function fstr_contact_init( contact, hecMESH, cparam, myrank )
157  type(tcontact), intent(inout) :: contact
158  type(hecmwst_local_mesh), pointer :: hecmesh
159  type(tcontactparam), target :: cparam
160  integer(kint),intent(in),optional :: myrank
161 
162  integer :: i, j, is, ie, cgrp, nsurf, nslave, ic, ic_type, iss, nn, ii
163  integer :: count,nodeid
164 
165  fstr_contact_init = .false.
166 
167  contact%cparam => cparam
168 
169  ! master surface
170  cgrp = contact%surf_id2
171  if( cgrp<=0 ) return
172  is= hecmesh%surf_group%grp_index(cgrp-1) + 1
173  ie= hecmesh%surf_group%grp_index(cgrp )
174 
175  if(present(myrank)) then
176  ! PARA_CONTACT
177  count = 0
178  do i=is,ie
179  ic = hecmesh%surf_group%grp_item(2*i-1)
180  if(hecmesh%elem_ID(ic*2) /= myrank) cycle
181  count = count + 1
182  enddo
183  allocate( contact%master(count) )
184  count = 0
185  do i=is,ie
186  ic = hecmesh%surf_group%grp_item(2*i-1)
187  if(hecmesh%elem_ID(ic*2) /= myrank) cycle
188  count = count + 1
189  nsurf = hecmesh%surf_group%grp_item(2*i)
190  ic_type = hecmesh%elem_type(ic)
191  call initialize_surf( ic, ic_type, nsurf, contact%master(count) )
192  iss = hecmesh%elem_node_index(ic-1)
193  do j=1, size( contact%master(count)%nodes )
194  nn = contact%master(count)%nodes(j)
195  contact%master(count)%nodes(j) = hecmesh%elem_node_item( iss+nn )
196  enddo
197  enddo
198 
199  else
200  ! not PARA_CONTACT
201  allocate( contact%master(ie-is+1) )
202  do i=is,ie
203  ic = hecmesh%surf_group%grp_item(2*i-1)
204  nsurf = hecmesh%surf_group%grp_item(2*i)
205  ic_type = hecmesh%elem_type(ic)
206  call initialize_surf( ic, ic_type, nsurf, contact%master(i-is+1) )
207  iss = hecmesh%elem_node_index(ic-1)
208  do j=1, size( contact%master(i-is+1)%nodes )
209  nn = contact%master(i-is+1)%nodes(j)
210  contact%master(i-is+1)%nodes(j) = hecmesh%elem_node_item( iss+nn )
211  enddo
212  enddo
213 
214  endif
215 
216  call update_surface_reflen( contact%master, hecmesh%node )
217 
218  cgrp = contact%surf_id1
219  if( cgrp<=0 ) return
220  is= hecmesh%node_group%grp_index(cgrp-1) + 1
221  ie= hecmesh%node_group%grp_index(cgrp )
222  nslave = 0
223  do i=is,ie
224  nodeid = hecmesh%global_node_ID(hecmesh%node_group%grp_item(i))
225  if(present(myrank)) then
226  ! PARA_CONTACT
227  nslave = nslave + 1
228  else
229  ! not PARA_CONTACT
230  if( hecmesh%node_group%grp_item(i) <= hecmesh%nn_internal) then
231  nslave = nslave + 1
232  endif
233  endif
234  enddo
235  allocate( contact%slave(nslave) )
236  ii = 0
237  do i=is,ie
238  if(.not.present(myrank)) then
239  ! not PARA_CONTACT
240  if( hecmesh%node_group%grp_item(i) > hecmesh%nn_internal) cycle
241  endif
242  ii = ii + 1
243  contact%slave(ii) = hecmesh%node_group%grp_item(i)
244  enddo
245 
246  ! contact state
247  allocate( contact%states(nslave) )
248  do i=1,nslave
249  call contact_state_init( contact%states(i) )
250  enddo
251 
252  ! neighborhood of surface group
253  call update_surface_box_info( contact%master, hecmesh%node )
254  call bucketdb_init( contact%master_bktDB )
255  call update_surface_bucket_info( contact%master, contact%master_bktDB )
256  call find_surface_neighbor( contact%master, contact%master_bktDB )
257 
258  ! initialize contact communication table
259  call hecmw_contact_comm_init( contact%comm, hecmesh, 1, nslave, contact%slave )
260 
261  contact%symmetric = .true.
262  fstr_contact_init = .true.
263  end function
264 
266  logical function fstr_embed_init( embed, hecMESH, cparam, myrank )
267  type(tcontact), intent(inout) :: embed
268  type(hecmwst_local_mesh), pointer :: hecmesh
269  type(tcontactparam), target :: cparam
270  integer(kint),intent(in),optional :: myrank
271 
272  integer :: i, j, is, ie, cgrp, nsurf, nslave, ic, ic_type, iss, nn, ii
273  integer :: count,nodeid
274 
275  fstr_embed_init = .false.
276 
277  embed%cparam => cparam
278 
279  ! master surface
280  cgrp = embed%surf_id2
281  if( cgrp<=0 ) return
282  is= hecmesh%elem_group%grp_index(cgrp-1) + 1
283  ie= hecmesh%elem_group%grp_index(cgrp )
284 
285  if(present(myrank)) then
286  ! PARA_CONTACT
287  count = 0
288  do i=is,ie
289  ic = hecmesh%elem_group%grp_item(i)
290  if(hecmesh%elem_ID(ic*2) /= myrank) cycle
291  count = count + 1
292  enddo
293  allocate( embed%master(count) )
294  count = 0
295  do i=is,ie
296  ic = hecmesh%elem_group%grp_item(i)
297  if(hecmesh%elem_ID(ic*2) /= myrank) cycle
298  count = count + 1
299  ic_type = hecmesh%elem_type(ic)
300  call initialize_surf( ic, ic_type, 0, embed%master(count) )
301  iss = hecmesh%elem_node_index(ic-1)
302  do j=1, size( embed%master(count)%nodes )
303  nn = embed%master(count)%nodes(j)
304  embed%master(count)%nodes(j) = hecmesh%elem_node_item( iss+nn )
305  enddo
306  enddo
307 
308  else
309  ! not PARA_CONTACT
310  allocate( embed%master(ie-is+1) )
311  do i=is,ie
312  ic = hecmesh%elem_group%grp_item(i)
313  ic_type = hecmesh%elem_type(ic)
314  call initialize_surf( ic, ic_type, 0, embed%master(i-is+1) )
315  iss = hecmesh%elem_node_index(ic-1)
316  do j=1, size( embed%master(i-is+1)%nodes )
317  nn = embed%master(i-is+1)%nodes(j)
318  embed%master(i-is+1)%nodes(j) = hecmesh%elem_node_item( iss+nn )
319  enddo
320  enddo
321 
322  endif
323 
324  ! slave surface
325  cgrp = embed%surf_id1
326  if( cgrp<=0 ) return
327  is= hecmesh%node_group%grp_index(cgrp-1) + 1
328  ie= hecmesh%node_group%grp_index(cgrp )
329  nslave = 0
330  do i=is,ie
331  nodeid = hecmesh%global_node_ID(hecmesh%node_group%grp_item(i))
332  if(present(myrank)) then
333  ! PARA_CONTACT
334  nslave = nslave + 1
335  else
336  ! not PARA_CONTACT
337  if( hecmesh%node_group%grp_item(i) <= hecmesh%nn_internal) then
338  nslave = nslave + 1
339  endif
340  endif
341  enddo
342  allocate( embed%slave(nslave) )
343  ii = 0
344  do i=is,ie
345  if(.not.present(myrank)) then
346  ! not PARA_CONTACT
347  if( hecmesh%node_group%grp_item(i) > hecmesh%nn_internal) cycle
348  endif
349  ii = ii + 1
350  embed%slave(ii) = hecmesh%node_group%grp_item(i)
351  enddo
352 
353  ! embed state
354  allocate( embed%states(nslave) )
355  do i=1,nslave
356  call contact_state_init( embed%states(i) )
357  enddo
358 
359  ! neighborhood of surface group
360  call update_surface_box_info( embed%master, hecmesh%node )
361  call bucketdb_init( embed%master_bktDB )
362  call update_surface_bucket_info( embed%master, embed%master_bktDB )
363  call find_surface_neighbor( embed%master, embed%master_bktDB )
364 
365  ! initialize contact communication table
366  call hecmw_contact_comm_init( embed%comm, hecmesh, 1, nslave, embed%slave )
367 
368  embed%symmetric = .true.
369  fstr_embed_init = .true.
370  end function
371 
372  function check_apply_contact_if( contact_if, contacts )
373  type(tcontactinterference), intent(inout) :: contact_if
374  type(tcontact) :: contacts(:)
375 
376  integer :: i, j
377  logical :: isfind
378  integer(kind=kint) :: check_apply_contact_if
379 
381  ! if contact pair exist?
382  isfind = .false.
383  do i = 1, size(contacts)
384  if( contacts(i)%pair_name == contact_if%cp_name ) then
385  contacts(i)%if_type = contact_if%if_type
386  contacts(i)%if_etime = contact_if%etime
387  contacts(i)%initial_pos = contact_if%initial_pos
388  contacts(i)%end_pos = contact_if%end_pos
389  do j = 1, size(contacts(i)%states)
390  contacts(i)%states(j)%interference_flag = contact_if%if_type
391  contacts(i)%states(j)%init_pos = contact_if%initial_pos
392  contacts(i)%states(j)%end_pos = contact_if%end_pos
393  if( contact_if%if_type /= c_if_slave )then
394  contacts(i)%states(j)%time_factor = (contact_if%end_pos - contact_if%initial_pos) / contact_if%etime
395  else
396  contacts(i)%states(j)%time_factor = contact_if%etime
397  end if
398  end do
399  isfind = .true.
400  check_apply_contact_if = 0; return
401  endif
402  enddo
403  if( .not. isfind ) return;
405 
406  end function
407 
409  subroutine clear_contact_state( contact )
410  type(tcontact), intent(inout) :: contact
411  integer :: i
412  if( .not. associated(contact%states) ) return
413  do i=1,size( contact%states )
414  contact%states(i)%state = -1
415  enddo
416  end subroutine
417 
419  logical function is_active_contact( acgrp, contact )
420  integer, intent(in) :: acgrp(:)
421  type(tcontact), intent(in) :: contact
422  if( any( acgrp==contact%group ) ) then
423  is_active_contact = .true.
424  else
425  is_active_contact = .false.
426  endif
427  end function
428 
432  subroutine scan_contact_state( flag_ctAlgo, contact, currpos, currdisp, ndforce, infoCTChange, &
433  nodeID, elemID, is_init, active, mu, B )
434  character(len=9), intent(in) :: flag_ctAlgo
435  type( tcontact ), intent(inout) :: contact
436  type( fstr_info_contactchange ), intent(inout) :: infoCTChange
437  real(kind=kreal), intent(in) :: currpos(:)
438  real(kind=kreal), intent(in) :: currdisp(:)
439  real(kind=kreal), intent(in) :: ndforce(:)
440  integer(kind=kint), intent(in) :: nodeID(:)
441  integer(kind=kint), intent(in) :: elemID(:)
442  logical, intent(in) :: is_init
443  logical, intent(out) :: active
444  real(kind=kreal), intent(in) :: mu
445  real(kind=kreal), optional, target :: b(:)
446 
447  real(kind=kreal) :: distclr
448  integer(kind=kint) :: slave, id, etype
449  integer(kind=kint) :: nn, i, j, iSS, nactive
450  real(kind=kreal) :: coord(3), elem(3, l_max_elem_node )
451  real(kind=kreal) :: nlforce, slforce(3)
452  logical :: isin
453  integer(kind=kint), allocatable :: contact_surf(:), states_prev(:)
454  !
455  integer, pointer :: indexMaster(:),indexCand(:)
456  integer :: nMaster,idm,nMasterMax,bktID,nCand
457  logical :: is_cand, is_present_B
458  real(kind=kreal), pointer :: bp(:)
459 
460  if( is_init ) then
461  distclr = contact%cparam%DISTCLR_INIT
462  else
463  distclr = contact%cparam%DISTCLR_FREE
464  if( contact%algtype == contacttied ) then
465  active = .false.
466  do i= 1, size(contact%slave)
467  if( contact%states(i)%state==contactstick ) then
468  active = .true.
469  exit
470  endif
471  enddo
472  endif
473  endif
474 
475  allocate(contact_surf(size(nodeid)))
476  allocate(states_prev(size(contact%slave)))
477  contact_surf(:) = huge(0)
478  do i = 1, size(contact%slave)
479  states_prev(i) = contact%states(i)%state
480  enddo
481 
482  call update_surface_box_info( contact%master, currpos )
483  call update_surface_bucket_info( contact%master, contact%master_bktDB )
484 
485  ! for gfortran-10: optional parameter seems not allowed within omp parallel
486  is_present_b = present(b)
487  if( is_present_b ) bp => b
488 
489  !$omp parallel do &
490  !$omp& default(none) &
491  !$omp& private(i,slave,slforce,id,nlforce,coord,indexMaster,nMaster,nn,j,iSS,elem,is_cand,idm,etype,isin, &
492  !$omp& bktID,nCand,indexCand) &
493  !$omp& firstprivate(nMasterMax,is_present_B) &
494  !$omp& shared(contact,ndforce,flag_ctAlgo,infoCTChange,currpos,currdisp,mu,nodeID,elemID,Bp,distclr,contact_surf,is_init) &
495  !$omp& reduction(.or.:active) &
496  !$omp& schedule(dynamic,1)
497  do i= 1, size(contact%slave)
498  if(contact%if_type /= 0) call set_shrink_factor(contact%ctime, contact%states(i), contact%if_etime, contact%if_type)
499  slave = contact%slave(i)
500  if( contact%states(i)%state==contactstick .or. contact%states(i)%state==contactslip ) then
501  slforce(1:3)=ndforce(3*slave-2:3*slave)
502  id = contact%states(i)%surface
503  nlforce = contact%states(i)%multiplier(1)
504 
505  ! update direction of TIED contact
506  if( contact%algtype == contacttied ) then
507  call update_direction( i, contact, currpos )
508  if (.not.is_init) cycle
509  endif
510 
511  if( nlforce < contact%cparam%TENSILE_FORCE ) then
512  contact%states(i)%state = contactfree
513  contact%states(i)%multiplier(:) = 0.d0
514  write(*,'(A,i10,A,i10,A,e12.3)') "Node",nodeid(slave)," free from contact with element", &
515  elemid(contact%master(id)%eid), " with tensile force ", nlforce
516  cycle
517  endif
518  if( contact%algtype /= contactfslid .or. (.not. is_present_b) ) then ! small slide problem
519  contact_surf(contact%slave(i)) = -elemid(contact%master(id)%eid)
520  else
521  call track_contact_position( flag_ctalgo, i, contact, currpos, currdisp, mu, infoctchange, nodeid, elemid, bp )
522  if( contact%states(i)%state /= contactfree ) then
523  id = contact%states(i)%surface
524  contact_surf(contact%slave(i)) = -elemid(contact%master(id)%eid)
525  endif
526  endif
527 
528  else if( contact%states(i)%state==contactfree ) then
529  if( contact%algtype == contacttied .and. .not. is_init ) cycle
530  coord(:) = currpos(3*slave-2:3*slave)
531 
532  ! get master candidates from bucketDB
533  bktid = bucketdb_getbucketid(contact%master_bktDB, coord)
534  ncand = bucketdb_getnumcand(contact%master_bktDB, bktid)
535  if (ncand == 0) cycle
536  allocate(indexcand(ncand))
537  call bucketdb_getcand(contact%master_bktDB, bktid, ncand, indexcand)
538 
539  nmastermax = ncand
540  allocate(indexmaster(nmastermax))
541  nmaster = 0
542 
543  ! narrow down candidates
544  do idm= 1, ncand
545  id = indexcand(idm)
546  if (.not. is_in_surface_box( contact%master(id), coord(1:3), contact%cparam%BOX_EXP_RATE )) cycle
547  nmaster = nmaster + 1
548  indexmaster(nmaster) = id
549  enddo
550  deallocate(indexcand)
551 
552  if(nmaster == 0) then
553  deallocate(indexmaster)
554  cycle
555  endif
556 
557  do idm = 1,nmaster
558  id = indexmaster(idm)
559  etype = contact%master(id)%etype
560  nn = size( contact%master(id)%nodes )
561  do j=1,nn
562  iss = contact%master(id)%nodes(j)
563  elem(1:3,j)=currpos(3*iss-2:3*iss)
564  enddo
565  call project_point2element( coord,etype,nn,elem,contact%master(id)%reflen,contact%states(i), &
566  isin,distclr,localclr=contact%cparam%CLEARANCE )
567  if( .not. isin ) cycle
568  contact%states(i)%surface = id
569  contact%states(i)%multiplier(:) = 0.d0
570  iss = isinsideelement( etype, contact%states(i)%lpos(1:2), contact%cparam%CLR_CAL_NORM )
571  if( iss>0 ) &
572  call cal_node_normal( id, iss, contact%master, currpos, contact%states(i)%lpos(1:2), &
573  contact%states(i)%direction(:) )
574  contact_surf(contact%slave(i)) = elemid(contact%master(id)%eid)
575  write(*,'(A,i10,A,i10,A,f7.3,A,2f7.3,A,3f7.3,A,i6)') "Node",nodeid(slave)," contact with element", &
576  elemid(contact%master(id)%eid), &
577  " with distance ", contact%states(i)%distance," at ",contact%states(i)%lpos(1:2), &
578  " along direction ", contact%states(i)%direction," rank=",hecmw_comm_get_rank()
579  exit
580  enddo
581  deallocate(indexmaster)
582  endif
583  enddo
584  !$omp end parallel do
585 
586  if( contact%algtype == contacttied .and. .not. is_init ) then
587  deallocate(contact_surf)
588  deallocate(states_prev)
589  return
590  endif
591 
592  call hecmw_contact_comm_allreduce_i(contact%comm, contact_surf, hecmw_min)
593  nactive = 0
594  do i = 1, size(contact%slave)
595  if (contact%states(i)%state /= contactfree) then ! any slave in contact
596  id = contact%states(i)%surface
597  if (abs(contact_surf(contact%slave(i))) /= elemid(contact%master(id)%eid)) then ! that is in contact with other surface
598  contact%states(i)%state = contactfree ! should be freed
599  write(*,'(A,i10,A,i10,A,i6,A,i6,A)') "Node",nodeid(contact%slave(i))," contact with element", &
600  & elemid(contact%master(id)%eid), " in rank",hecmw_comm_get_rank()," freed due to duplication"
601  else
602  nactive = nactive + 1
603  endif
604  endif
605  if (states_prev(i) == contactfree .and. contact%states(i)%state /= contactfree) then
606  infoctchange%free2contact = infoctchange%free2contact + 1
607  elseif (states_prev(i) /= contactfree .and. contact%states(i)%state == contactfree) then
608  infoctchange%contact2free = infoctchange%contact2free + 1
609  endif
610  enddo
611  active = (nactive > 0)
612  deallocate(contact_surf)
613  deallocate(states_prev)
614  end subroutine scan_contact_state
615 
619  subroutine scan_embed_state( flag_ctAlgo, embed, currpos, currdisp, ndforce, infoCTChange, &
620  nodeID, elemID, is_init, active, mu, B )
621  character(len=9), intent(in) :: flag_ctAlgo
622  type( tcontact ), intent(inout) :: embed
623  type( fstr_info_contactchange ), intent(inout) :: infoCTChange
624  real(kind=kreal), intent(in) :: currpos(:)
625  real(kind=kreal), intent(in) :: currdisp(:)
626  real(kind=kreal), intent(in) :: ndforce(:)
627  integer(kind=kint), intent(in) :: nodeID(:)
628  integer(kind=kint), intent(in) :: elemID(:)
629  logical, intent(in) :: is_init
630  logical, intent(out) :: active
631  real(kind=kreal), intent(in) :: mu
632  real(kind=kreal), optional, target :: b(:)
633 
634  real(kind=kreal) :: distclr
635  integer(kind=kint) :: slave, id, etype
636  integer(kind=kint) :: nn, i, j, iSS, nactive
637  real(kind=kreal) :: coord(3), elem(3, l_max_elem_node )
638  real(kind=kreal) :: nlforce, slforce(3)
639  logical :: isin
640  integer(kind=kint), allocatable :: contact_surf(:), states_prev(:)
641  !
642  integer, pointer :: indexMaster(:),indexCand(:)
643  integer :: nMaster,idm,nMasterMax,bktID,nCand
644  logical :: is_cand, is_present_B
645  real(kind=kreal), pointer :: bp(:)
646 
647  if( is_init ) then
648  distclr = embed%cparam%DISTCLR_INIT
649  else
650  distclr = embed%cparam%DISTCLR_FREE
651  active = .false.
652  do i= 1, size(embed%slave)
653  if( embed%states(i)%state==contactstick ) then
654  active = .true.
655  exit
656  endif
657  enddo
658  return
659  endif
660 
661  allocate(contact_surf(size(nodeid)))
662  allocate(states_prev(size(embed%slave)))
663  contact_surf(:) = huge(0)
664  do i = 1, size(embed%slave)
665  states_prev(i) = embed%states(i)%state
666  enddo
667 
668  call update_surface_box_info( embed%master, currpos )
669  call update_surface_bucket_info( embed%master, embed%master_bktDB )
670 
671  ! for gfortran-10: optional parameter seems not allowed within omp parallel
672  is_present_b = present(b)
673  if( is_present_b ) bp => b
674 
675  !$omp parallel do &
676  !$omp& default(none) &
677  !$omp& private(i,slave,slforce,id,nlforce,coord,indexMaster,nMaster,nn,j,iSS,elem,is_cand,idm,etype,isin, &
678  !$omp& bktID,nCand,indexCand) &
679  !$omp& firstprivate(nMasterMax,is_present_B) &
680  !$omp& shared(embed,ndforce,flag_ctAlgo,infoCTChange,currpos,currdisp,mu,nodeID,elemID,Bp,distclr,contact_surf) &
681  !$omp& reduction(.or.:active) &
682  !$omp& schedule(dynamic,1)
683  do i= 1, size(embed%slave)
684  slave = embed%slave(i)
685  if( embed%states(i)%state==contactfree ) then
686  coord(:) = currpos(3*slave-2:3*slave)
687 
688  ! get master candidates from bucketDB
689  bktid = bucketdb_getbucketid(embed%master_bktDB, coord)
690  ncand = bucketdb_getnumcand(embed%master_bktDB, bktid)
691  if (ncand == 0) cycle
692  allocate(indexcand(ncand))
693  call bucketdb_getcand(embed%master_bktDB, bktid, ncand, indexcand)
694 
695  nmastermax = ncand
696  allocate(indexmaster(nmastermax))
697  nmaster = 0
698 
699  ! narrow down candidates
700  do idm= 1, ncand
701  id = indexcand(idm)
702  if (.not. is_in_surface_box( embed%master(id), coord(1:3), embed%cparam%BOX_EXP_RATE )) cycle
703  nmaster = nmaster + 1
704  indexmaster(nmaster) = id
705  enddo
706  deallocate(indexcand)
707 
708  if(nmaster == 0) then
709  deallocate(indexmaster)
710  cycle
711  endif
712 
713  do idm = 1,nmaster
714  id = indexmaster(idm)
715  etype = embed%master(id)%etype
716  if( mod(etype,10) == 2 ) etype = etype - 1 !search by 1st-order shape function
717  nn = getnumberofnodes(etype)
718  do j=1,nn
719  iss = embed%master(id)%nodes(j)
720  elem(1:3,j)=currpos(3*iss-2:3*iss)
721  enddo
722  call project_point2solidelement( coord,etype,nn,elem,embed%master(id)%reflen,embed%states(i), &
723  isin,distclr,localclr=embed%cparam%CLEARANCE )
724  if( .not. isin ) cycle
725  embed%states(i)%surface = id
726  embed%states(i)%multiplier(:) = 0.d0
727  contact_surf(embed%slave(i)) = elemid(embed%master(id)%eid)
728  write(*,'(A,i10,A,i10,A,3f7.3,A,i6)') "Node",nodeid(slave)," embeded to element", &
729  elemid(embed%master(id)%eid), " at ",embed%states(i)%lpos(:)," rank=",hecmw_comm_get_rank()
730  exit
731  enddo
732  deallocate(indexmaster)
733  endif
734  enddo
735  !$omp end parallel do
736 
737  call hecmw_contact_comm_allreduce_i(embed%comm, contact_surf, hecmw_min)
738  nactive = 0
739  do i = 1, size(embed%slave)
740  if (embed%states(i)%state /= contactfree) then ! any slave in contact
741  id = embed%states(i)%surface
742  if (abs(contact_surf(embed%slave(i))) /= elemid(embed%master(id)%eid)) then ! that is in contact with other surface
743  embed%states(i)%state = contactfree ! should be freed
744  write(*,'(A,i10,A,i10,A,i6,A,i6,A)') "Node",nodeid(embed%slave(i))," contact with element", &
745  & elemid(embed%master(id)%eid), " in rank",hecmw_comm_get_rank()," freed due to duplication"
746  else
747  nactive = nactive + 1
748  endif
749  endif
750  if (states_prev(i) == contactfree .and. embed%states(i)%state /= contactfree) then
751  infoctchange%free2contact = infoctchange%free2contact + 1
752  elseif (states_prev(i) /= contactfree .and. embed%states(i)%state == contactfree) then
753  infoctchange%contact2free = infoctchange%contact2free + 1
754  endif
755  enddo
756  active = (nactive > 0)
757  deallocate(contact_surf)
758  deallocate(states_prev)
759  end subroutine scan_embed_state
760 
761 
763  subroutine cal_node_normal( csurf, isin, surf, currpos, lpos, normal )
765  integer, intent(in) :: csurf
766  integer, intent(in) :: isin
767  type(tsurfelement), intent(in) :: surf(:)
768  real(kind=kreal), intent(in) :: currpos(:)
769  real(kind=kreal), intent(in) :: lpos(:)
770  real(kind=kreal), intent(out) :: normal(3)
771  integer(kind=kint) :: cnode, i, j, cnt, nd1, gn, etype, iSS, nn,cgn
772  real(kind=kreal) :: cnpos(2), elem(3, l_max_elem_node )
773  integer(kind=kint) :: cnode1, cnode2, gn1, gn2, nsurf, cgn1, cgn2, isin_n
774  real(kind=kreal) :: x=0, normal_n(3), lpos_n(2)
775 
776  if( 1 <= isin .and. isin <= 4 ) then ! corner
777  cnode = isin
778  gn = surf(csurf)%nodes(cnode)
779  etype = surf(csurf)%etype
780  call getvertexcoord( etype, cnode, cnpos )
781  nn = size( surf(csurf)%nodes )
782  do j=1,nn
783  iss = surf(csurf)%nodes(j)
784  elem(1:3,j)=currpos(3*iss-2:3*iss)
785  enddo
786  normal = surfacenormal( etype, nn, cnpos, elem )
787  cnt = 1
788  do i=1,surf(csurf)%n_neighbor
789  nd1 = surf(csurf)%neighbor(i)
790  nn = size( surf(nd1)%nodes )
791  etype = surf(nd1)%etype
792  cgn = 0
793  do j=1,nn
794  iss = surf(nd1)%nodes(j)
795  elem(1:3,j)=currpos(3*iss-2:3*iss)
796  if( iss==gn ) cgn=j
797  enddo
798  if( cgn>0 ) then
799  call getvertexcoord( etype, cgn, cnpos )
800  !normal = normal+SurfaceNormal( etype, nn, cnpos, elem )
801  normal_n = surfacenormal( etype, nn, cnpos, elem )
802  normal = normal+normal_n
803  cnt = cnt+1
804  endif
805  enddo
806  !normal = normal/cnt !!-???
807  elseif( 12 <= isin .and. isin <= 41 ) then ! edge
808  cnode1 = isin / 10
809  cnode2 = mod(isin, 10)
810  gn1 = surf(csurf)%nodes(cnode1)
811  gn2 = surf(csurf)%nodes(cnode2)
812  etype = surf(csurf)%etype
813  nn = size( surf(csurf)%nodes )
814  do j=1,nn
815  iss = surf(csurf)%nodes(j)
816  elem(1:3,j)=currpos(3*iss-2:3*iss)
817  enddo
818  normal = surfacenormal( etype, nn, lpos, elem )
819  select case (etype)
820  case (fe_tri3n, fe_tri6n, fe_tri6nc)
821  if ( isin==12 ) then; x=lpos(2)-lpos(1)
822  elseif( isin==23 ) then; x=1.d0-2.d0*lpos(2)
823  elseif( isin==31 ) then; x=2.d0*lpos(1)-1.d0
824  else; stop "Error: cal_node_normal: invalid isin"
825  endif
826  case (fe_quad4n, fe_quad8n)
827  if ( isin==12 ) then; x=lpos(1)
828  elseif( isin==23 ) then; x=lpos(2)
829  elseif( isin==34 ) then; x=-lpos(1)
830  elseif( isin==41 ) then; x=-lpos(2)
831  else; stop "Error: cal_node_normal: invalid isin"
832  endif
833  end select
834  ! find neighbor surf that includes cnode1 and cnode2
835  nsurf = 0
836  neib_loop: do i=1, surf(csurf)%n_neighbor
837  nd1 = surf(csurf)%neighbor(i)
838  nn = size( surf(nd1)%nodes )
839  etype = surf(nd1)%etype
840  cgn1 = 0
841  cgn2 = 0
842  do j=1,nn
843  iss = surf(nd1)%nodes(j)
844  elem(1:3,j)=currpos(3*iss-2:3*iss)
845  if( iss==gn1 ) cgn1=j
846  if( iss==gn2 ) cgn2=j
847  enddo
848  if( cgn1>0 .and. cgn2>0 ) then
849  nsurf = nd1
850  isin_n = 10*cgn2 + cgn1
851  x = -x
852  select case (etype)
853  case (fe_tri3n, fe_tri6n, fe_tri6nc)
854  if ( isin_n==12 ) then; lpos_n(1)=0.5d0*(1.d0-x); lpos_n(2)=0.5d0*(1.d0+x)
855  elseif( isin_n==23 ) then; lpos_n(1)=0.d0; lpos_n(2)=0.5d0*(1.d0-x)
856  elseif( isin_n==31 ) then; lpos_n(1)=0.5d0*(1.d0+x); lpos_n(2)=0.d0
857  else; stop "Error: cal_node_normal: invalid isin_n"
858  endif
859  case (fe_quad4n, fe_quad8n)
860  if ( isin_n==12 ) then; lpos_n(1)= x; lpos_n(2)=-1.d0
861  elseif( isin_n==23 ) then; lpos_n(1)= 1.d0; lpos_n(2)= x
862  elseif( isin_n==34 ) then; lpos_n(1)=-x; lpos_n(2)= 1.d0
863  elseif( isin_n==41 ) then; lpos_n(1)=-1.d0; lpos_n(2)=-x
864  else; stop "Error: cal_node_normal: invalid isin_n"
865  endif
866  end select
867  !normal = normal + SurfaceNormal( etype, nn, lpos_n, elem )
868  normal_n = surfacenormal( etype, nn, lpos_n, elem )
869  normal = normal+normal_n
870  exit neib_loop
871  endif
872  enddo neib_loop
873  !if( nsurf==0 ) write(0,*) "Warning: cal_node_normal: neighbor surf not found"
874  !normal = normal/2
875  endif
876  normal = normal/ dsqrt( dot_product( normal, normal ) )
877  end subroutine cal_node_normal
878 
880  subroutine track_contact_position( flag_ctAlgo, nslave, contact, currpos, currdisp, mu, infoCTChange, nodeID, elemID, B )
881  character(len=9), intent(in) :: flag_ctAlgo
882  integer, intent(in) :: nslave
883  type( tcontact ), intent(inout) :: contact
884  type( fstr_info_contactchange ), intent(inout) :: infoCTChange
885  real(kind=kreal), intent(in) :: currpos(:)
886  real(kind=kreal), intent(in) :: currdisp(:)
887  real(kind=kreal), intent(in) :: mu
888  integer(kind=kint), intent(in) :: nodeID(:)
889  integer(kind=kint), intent(in) :: elemID(:)
890  real(kind=kreal), intent(inout) :: b(:)
891 
892  integer(kind=kint) :: slave, sid0, sid, etype
893  integer(kind=kint) :: nn, i, j, iSS
894  real(kind=kreal) :: coord(3), elem(3, l_max_elem_node ), elem0(3, l_max_elem_node )
895  logical :: isin
896  real(kind=kreal) :: opos(2), odirec(3)
897  integer(kind=kint) :: bktID, nCand, idm
898  integer(kind=kint), allocatable :: indexCand(:)
899 
900  sid = 0
901 
902  slave = contact%slave(nslave)
903  coord(:) = currpos(3*slave-2:3*slave)
905  sid0 = contact%states(nslave)%surface
906  opos = contact%states(nslave)%lpos(1:2)
907  odirec = contact%states(nslave)%direction
908  etype = contact%master(sid0)%etype
909  nn = getnumberofnodes( etype )
910  do j=1,nn
911  iss = contact%master(sid0)%nodes(j)
912  elem(1:3,j)=currpos(3*iss-2:3*iss)
913  elem0(1:3,j)=currpos(3*iss-2:3*iss)-currdisp(3*iss-2:3*iss)
914  enddo
915  call project_point2element( coord,etype,nn,elem,contact%master(sid0)%reflen,contact%states(nslave), &
916  isin,contact%cparam%DISTCLR_NOCHECK,contact%states(nslave)%lpos(1:2),contact%cparam%CLR_SAME_ELEM )
917  if( .not. isin ) then
918  do i=1, contact%master(sid0)%n_neighbor
919  sid = contact%master(sid0)%neighbor(i)
920  etype = contact%master(sid)%etype
921  nn = getnumberofnodes( etype )
922  do j=1,nn
923  iss = contact%master(sid)%nodes(j)
924  elem(1:3,j)=currpos(3*iss-2:3*iss)
925  enddo
926  call project_point2element( coord,etype,nn,elem,contact%master(sid)%reflen,contact%states(nslave), &
927  isin,contact%cparam%DISTCLR_NOCHECK,localclr=contact%cparam%CLEARANCE )
928  if( isin ) then
929  contact%states(nslave)%surface = sid
930  exit
931  endif
932  enddo
933  endif
934 
935  if( .not. isin ) then ! such case is considered to rarely or never occur
936  write(*,*) 'Warning: contact moved beyond neighbor elements'
937  ! get master candidates from bucketDB
938  bktid = bucketdb_getbucketid(contact%master_bktDB, coord)
939  ncand = bucketdb_getnumcand(contact%master_bktDB, bktid)
940  if (ncand > 0) then
941  allocate(indexcand(ncand))
942  call bucketdb_getcand(contact%master_bktDB, bktid, ncand, indexcand)
943  do idm= 1, ncand
944  sid = indexcand(idm)
945  if( sid==sid0 ) cycle
946  if( associated(contact%master(sid0)%neighbor) ) then
947  if( any(sid==contact%master(sid0)%neighbor(:)) ) cycle
948  endif
949  if (.not. is_in_surface_box( contact%master(sid), coord(1:3), contact%cparam%BOX_EXP_RATE )) cycle
950  etype = contact%master(sid)%etype
951  nn = size( contact%master(sid)%nodes )
952  do j=1,nn
953  iss = contact%master(sid)%nodes(j)
954  elem(1:3,j)=currpos(3*iss-2:3*iss)
955  enddo
956  call project_point2element( coord,etype,nn,elem,contact%master(sid)%reflen,contact%states(nslave), &
957  isin,contact%cparam%DISTCLR_NOCHECK,localclr=contact%cparam%CLEARANCE )
958  if( isin ) then
959  contact%states(nslave)%surface = sid
960  exit
961  endif
962  enddo
963  deallocate(indexcand)
964  endif
965  endif
966 
967  if( isin ) then
968  if( contact%states(nslave)%surface==sid0 ) then
969  if(any(dabs(contact%states(nslave)%lpos(1:2)-opos(:)) >= contact%cparam%CLR_DIFFLPOS)) then
970  !$omp atomic
971  infoctchange%contact2difflpos = infoctchange%contact2difflpos + 1
972  endif
973  else
974  write(*,'(A,i10,A,i10,A,f7.3,A,2f7.3)') "Node",nodeid(slave)," move to contact with", &
975  elemid(contact%master(sid)%eid), " with distance ", &
976  contact%states(nslave)%distance," at ",contact%states(nslave)%lpos(1:2)
977  !$omp atomic
978  infoctchange%contact2neighbor = infoctchange%contact2neighbor + 1
979  if( flag_ctalgo=='ALagrange' ) &
980  call reset_contact_force( contact, currpos, nslave, sid0, opos, odirec, b )
981  endif
982  if( flag_ctalgo=='SLagrange' ) call update_tangentforce(etype,nn,elem0,elem,contact%states(nslave))
983  iss = isinsideelement( etype, contact%states(nslave)%lpos(1:2), contact%cparam%CLR_CAL_NORM )
984  if( iss>0 ) &
985  call cal_node_normal( contact%states(nslave)%surface, iss, contact%master, currpos, &
986  contact%states(nslave)%lpos(1:2), contact%states(nslave)%direction(:) )
987  else if( .not. isin ) then
988  write(*,'(A,i10,A)') "Node",nodeid(slave)," move out of contact"
989  contact%states(nslave)%state = contactfree
990  contact%states(nslave)%multiplier(:) = 0.d0
991  endif
992 
993  end subroutine track_contact_position
994 
996  subroutine update_direction( nslave, contact, currpos )
997  integer, intent(in) :: nslave
998  type( tcontact ), intent(inout) :: contact
999  real(kind=kreal), intent(in) :: currpos(:)
1000 
1001  integer(kind=kint) :: slave, sid0, sid, etype
1002  integer(kind=kint) :: nn, i, j, iSS
1003  real(kind=kreal) :: coord(3), elem(3, l_max_elem_node ), elem0(3, l_max_elem_node )
1004  logical :: isin
1005  real(kind=kreal) :: opos(2), odirec(3)
1006  integer(kind=kint) :: bktID, nCand, idm
1007  integer(kind=kint), allocatable :: indexCand(:)
1008  type(tcontactstate) :: cstate_tmp
1009 
1010  sid = 0
1011 
1012  slave = contact%slave(nslave)
1013  coord(:) = currpos(3*slave-2:3*slave)
1015  sid0 = contact%states(nslave)%surface
1016  opos = contact%states(nslave)%lpos(1:2)
1017  odirec = contact%states(nslave)%direction
1018  etype = contact%master(sid0)%etype
1019  nn = getnumberofnodes( etype )
1020  do j=1,nn
1021  iss = contact%master(sid0)%nodes(j)
1022  elem(1:3,j)=currpos(3*iss-2:3*iss)
1023  enddo
1024 
1025  cstate_tmp%state = contact%states(nslave)%state
1026  cstate_tmp%surface = sid0
1027  call project_point2element( coord,etype,nn,elem,contact%master(sid0)%reflen,cstate_tmp, &
1028  isin,contact%cparam%DISTCLR_NOCHECK,contact%states(nslave)%lpos,contact%cparam%CLR_SAME_ELEM )
1029 
1030  if( isin ) then
1031  iss = isinsideelement( etype, cstate_tmp%lpos, contact%cparam%CLR_CAL_NORM )
1032  if( iss>0 ) &
1033  call cal_node_normal( cstate_tmp%surface, iss, contact%master, currpos, &
1034  cstate_tmp%lpos, cstate_tmp%direction(:) )
1035  endif
1036 
1037  contact%states(nslave)%direction = cstate_tmp%direction
1038 
1039  end subroutine
1040 
1041 
1044  subroutine reset_contact_force( contact, currpos, lslave, omaster, opos, odirec, B )
1045  type( tcontact ), intent(inout) :: contact
1046  real(kind=kreal), intent(in) :: currpos(:)
1047  integer, intent(in) :: lslave
1048  integer, intent(in) :: omaster
1049  real(kind=kreal), intent(in) :: opos(2)
1050  real(kind=kreal), intent(in) :: odirec(3)
1051  real(kind=kreal), intent(inout) :: b(:)
1052 
1053  integer(kind=kint) :: slave, etype, master
1054  integer(kind=kint) :: nn, j, iSS
1055  real(kind=kreal) :: nrlforce, fcoeff, tangent(3,2)
1056  real(kind=kreal) :: elemcrd(3, l_max_elem_node )
1057  real(kind=kreal) :: shapefunc(l_max_surface_node)
1058  real(kind=kreal) :: metric(2,2), dispmat(2,l_max_elem_node*3+3)
1059  real(kind=kreal) :: fric(2), f3(l_max_elem_node*3+3)
1060  integer(kind=kint) :: i, idx0
1061 
1062  slave = contact%slave(lslave)
1063  fcoeff = contact%fcoeff
1064  ! clear contact force in former contacted element
1065  nrlforce = contact%states(lslave)%multiplier(1)
1066  b(3*slave-2:3*slave) = b(3*slave-2:3*slave)+nrlforce*odirec
1067  nn = size( contact%master(omaster)%nodes )
1068  etype = contact%master(omaster)%etype
1069  call getshapefunc( etype, opos(:), shapefunc )
1070  do j=1,nn
1071  iss = contact%master(omaster)%nodes(j)
1072  ! B(3*iSS-2:3*iSS) = B(3*iSS-2:3*iSS)-nrlforce*shapefunc(j)*odirec
1073  idx0 = 3*(iss-1)
1074  do i=1,3
1075  !$omp atomic
1076  b(idx0+i) = b(idx0+i)-nrlforce*shapefunc(j)*odirec(i)
1077  enddo
1078  enddo
1079  if( fcoeff/=0.d0 ) then
1080  do j=1,nn
1081  iss = contact%master(omaster)%nodes(j)
1082  elemcrd(:,j) = currpos(3*iss-2:3*iss)
1083  enddo
1084  call dispincrematrix( opos(:), etype, nn, elemcrd, tangent, &
1085  metric, dispmat )
1086  fric(1:2) = contact%states(lslave)%multiplier(2:3)
1087  f3(:) = fric(1)*dispmat(1,:)+fric(2)*dispmat(2,:)
1088  b(3*slave-2:3*slave) = b(3*slave-2:3*slave)+f3(1:3)
1089  do j=1,nn
1090  iss = contact%master(omaster)%nodes(j)
1091  ! B(3*iSS-2:3*iSS) = B(3*iSS-2:3*iSS)+f3(3*j+1:3*j+3)
1092  idx0 = 3*(iss-1)
1093  do i=1,3
1094  !$omp atomic
1095  b(idx0+i) = b(idx0+i)+f3(3*j+i)
1096  enddo
1097  enddo
1098  endif
1099 
1100  ! reset contact force in new contacting element
1101  master = contact%states(lslave)%surface
1102  nn = size( contact%master(master)%nodes )
1103  etype = contact%master(master)%etype
1104  call getshapefunc( etype, contact%states(lslave)%lpos(1:2), shapefunc )
1105  b(3*slave-2:3*slave) = b(3*slave-2:3*slave)-nrlforce*contact%states(lslave)%direction
1106  do j=1,nn
1107  iss = contact%master(master)%nodes(j)
1108  ! B(3*iSS-2:3*iSS) = B(3*iSS-2:3*iSS)+nrlforce* &
1109  ! shapefunc(j)*contact%states(lslave)%direction
1110  idx0 = 3*(iss-1)
1111  do i=1,3
1112  !$omp atomic
1113  b(idx0+i) = b(idx0+i)+nrlforce* &
1114  shapefunc(j)*contact%states(lslave)%direction(i)
1115  enddo
1116  enddo
1117  if( fcoeff/=0.d0 ) then
1118  do j=1,nn
1119  iss = contact%master(master)%nodes(j)
1120  elemcrd(:,j) = currpos(3*iss-2:3*iss)
1121  enddo
1122  call dispincrematrix( contact%states(lslave)%lpos(1:2), etype, nn, elemcrd, tangent, &
1123  metric, dispmat )
1124  fric(1:2) = contact%states(lslave)%multiplier(2:3)
1125  f3(:) = fric(1)*dispmat(1,:)+fric(2)*dispmat(2,:)
1126  b(3*slave-2:3*slave) = b(3*slave-2:3*slave)-f3(1:3)
1127  do j=1,nn
1128  iss = contact%master(master)%nodes(j)
1129  ! B(3*iSS-2:3*iSS) = B(3*iSS-2:3*iSS)-f3(3*j+1:3*j+3)
1130  idx0 = 3*(iss-1)
1131  do i=1,3
1132  !$omp atomic
1133  b(idx0+i) = b(idx0+i)-f3(3*j+i)
1134  enddo
1135  enddo
1136  endif
1137 
1138  end subroutine reset_contact_force
1139 
1143  subroutine calcu_contact_force0( contact, coord, disp, ddisp, fcoeff, mu, &
1144  mut, B )
1145  type( tcontact ), intent(inout) :: contact
1146  real(kind=kreal), intent(in) :: coord(:)
1147  real(kind=kreal), intent(in) :: disp(:)
1148  real(kind=kreal), intent(in) :: ddisp(:)
1149  real(kind=kreal), intent(in) :: fcoeff
1150  real(kind=kreal), intent(in) :: mu, mut
1151  real(kind=kreal), intent(inout) :: b(:)
1152 
1153  integer(kind=kint) :: slave, etype, master
1154  integer(kind=kint) :: nn, i, j, iSS
1155  real(kind=kreal) :: nrlforce, elemdisp(3,l_max_elem_node), tangent(3,2)
1156  real(kind=kreal) :: dg(3), elemg(3), elemcrd(3, l_max_elem_node )
1157  real(kind=kreal) :: dgn, dxi(2), dxy(2), shapefunc(l_max_surface_node)
1158  real(kind=kreal) :: metric(2,2), dispmat(2,l_max_elem_node*3+3)
1159  real(kind=kreal) :: fric(2), f3(3*l_max_elem_node+3), edisp(3*l_max_elem_node+3)
1160 
1161  do i= 1, size(contact%slave)
1162  if( contact%states(i)%state==contactfree ) cycle ! not in contact
1163  slave = contact%slave(i)
1164  edisp(1:3) = ddisp(3*slave-2:3*slave)
1165  master = contact%states(i)%surface
1166 
1167  nn = size( contact%master(master)%nodes )
1168  etype = contact%master(master)%etype
1169  do j=1,nn
1170  iss = contact%master(master)%nodes(j)
1171  elemdisp(:,j) = ddisp(3*iss-2:3*iss)
1172  edisp(3*j+1:3*j+3) = ddisp(3*iss-2:3*iss)
1173  elemcrd(:,j) = coord(3*iss-2:3*iss)+disp(3*iss-2:3*iss)
1174  enddo
1175  call getshapefunc( etype, contact%states(i)%lpos(1:2), shapefunc )
1176 
1177  ! normal component
1178  elemg = 0.d0
1179  do j=1,nn
1180  elemg(:) = elemg(:)+shapefunc(j)*elemdisp(:,j)
1181  enddo
1182  dg(:) = ddisp(3*slave-2:3*slave) - elemg(:)
1183  dgn = dot_product( contact%states(i)%direction, dg )
1184  nrlforce = contact%states(i)%multiplier(1)- mu*(contact%states(i)%wkdist-dgn)
1185  b(3*slave-2:3*slave) = b(3*slave-2:3*slave)-nrlforce*contact%states(i)%direction
1186 
1187  do j=1,nn
1188  iss = contact%master(master)%nodes(j)
1189  b(3*iss-2:3*iss) = b(3*iss-2:3*iss)+nrlforce* &
1190  shapefunc(j)*contact%states(i)%direction
1191  enddo
1192 
1193  if( fcoeff==0.d0 ) cycle
1194  ! tangent component
1195  call dispincrematrix( contact%states(i)%lpos(1:2), etype, nn, elemcrd, tangent, &
1196  metric, dispmat )
1197  dxi(1) = dot_product( dispmat(1,1:nn*3+3), edisp(1:nn*3+3) )
1198  dxi(2) = dot_product( dispmat(2,1:nn*3+3), edisp(1:nn*3+3) )
1199  dxy(:) = matmul( metric, dxi )
1200  fric(1:2) = contact%states(i)%multiplier(2:3) + mut*dxy(1:2)
1201  f3(:) = fric(1)*dispmat(1,:)+fric(2)*dispmat(2,:)
1202 
1203  if( contact%states(i)%state==contactslip ) then
1204  dgn = dsqrt( f3(1)*f3(1)+f3(2)*f3(2)+f3(3)*f3(3) )
1205  f3(:) = f3(:)*fcoeff*contact%states(i)%multiplier(1)/dgn
1206  endif
1207  b(3*slave-2:3*slave) = b(3*slave-2:3*slave)-f3(1:3)
1208  do j=1,nn
1209  iss = contact%master(master)%nodes(j)
1210  b(3*iss-2:3*iss) = b(3*iss-2:3*iss)-f3(3*j+1:3*j+3)
1211  enddo
1212  enddo
1213  end subroutine calcu_contact_force0
1214 
1218  subroutine calcu_tied_force0( contact, disp, ddisp, mu, B )
1219  type( tcontact ), intent(inout) :: contact
1220  real(kind=kreal), intent(in) :: disp(:)
1221  real(kind=kreal), intent(in) :: ddisp(:)
1222  real(kind=kreal), intent(in) :: mu
1223  real(kind=kreal), intent(inout) :: b(:)
1224 
1225  integer(kind=kint) :: slave, etype, master
1226  integer(kind=kint) :: nn, i, j, iss
1227  real(kind=kreal) :: nrlforce(3)
1228  real(kind=kreal) :: dg(3)
1229  real(kind=kreal) :: shapefunc(l_max_surface_node)
1230  real(kind=kreal) :: edisp(3*l_max_elem_node+3)
1231 
1232  do i= 1, size(contact%slave)
1233  if( contact%states(i)%state==contactfree ) cycle ! not in contact
1234  slave = contact%slave(i)
1235  edisp(1:3) = disp(3*slave-2:3*slave)+ddisp(3*slave-2:3*slave)
1236  master = contact%states(i)%surface
1237 
1238  nn = size( contact%master(master)%nodes )
1239  etype = contact%master(master)%etype
1240  do j=1,nn
1241  iss = contact%master(master)%nodes(j)
1242  edisp(3*j+1:3*j+3) = disp(3*iss-2:3*iss)+ddisp(3*iss-2:3*iss)
1243  enddo
1244  call getshapefunc( etype, contact%states(i)%lpos(1:2), shapefunc )
1245 
1246  ! normal component
1247  dg(1:3) = edisp(1:3)
1248  do j=1,nn
1249  dg(1:3) = dg(1:3)-shapefunc(j)*edisp(3*j+1:3*j+3)
1250  enddo
1251 
1252  nrlforce(1:3) = contact%states(i)%multiplier(1:3)+mu*dg(1:3)
1253 
1254  b(3*slave-2:3*slave) = b(3*slave-2:3*slave)-nrlforce(1:3)
1255  do j=1,nn
1256  iss = contact%master(master)%nodes(j)
1257  b(3*iss-2:3*iss) = b(3*iss-2:3*iss) + shapefunc(j)*nrlforce(1:3)
1258  enddo
1259  enddo
1260 
1261  end subroutine
1262 
1265  subroutine update_contact_multiplier( contact, coord, disp, ddisp, fcoeff, mu, mut, &
1266  gnt, ctchanged )
1267  type( tcontact ), intent(inout) :: contact
1268  real(kind=kreal), intent(in) :: coord(:)
1269  real(kind=kreal), intent(in) :: disp(:)
1270  real(kind=kreal), intent(in) :: ddisp(:)
1271  real(kind=kreal), intent(in) :: fcoeff
1272  real(kind=kreal), intent(in) :: mu, mut
1273  real(kind=kreal), intent(out) :: gnt(2)
1274  logical, intent(inout) :: ctchanged
1275 
1276  integer(kind=kint) :: slave, etype, master
1277  integer(kind=kint) :: nn, i, j, iss, cnt
1278  real(kind=kreal) :: elemdisp(3,l_max_elem_node), tangent(3,2)
1279  real(kind=kreal) :: dg(3), elemg(3), elemcrd(3, l_max_elem_node )
1280  real(kind=kreal) :: dgn, dxi(2), dxy(2), shapefunc(l_max_surface_node)
1281  real(kind=kreal) :: metric(2,2), dispmat(2,l_max_elem_node*3+3)
1282  real(kind=kreal) :: lgnt(2), fric(2), f3(3*l_max_elem_node+3), edisp(3*l_max_elem_node+3)
1283 
1284  cnt =0; lgnt(:)=0.d0
1285  do i= 1, size(contact%slave)
1286  if( contact%states(i)%state==contactfree ) cycle ! not in contact
1287  slave = contact%slave(i)
1288  edisp(1:3) = ddisp(3*slave-2:3*slave)
1289  master = contact%states(i)%surface
1290 
1291  nn = size( contact%master(master)%nodes )
1292  etype = contact%master(master)%etype
1293  do j=1,nn
1294  iss = contact%master(master)%nodes(j)
1295  elemdisp(:,j) = ddisp(3*iss-2:3*iss)
1296  edisp(3*j+1:3*j+3) = ddisp(3*iss-2:3*iss)
1297  elemcrd(:,j) = coord(3*iss-2:3*iss)+disp(3*iss-2:3*iss)
1298  enddo
1299  call getshapefunc( etype, contact%states(i)%lpos(1:2), shapefunc )
1300 
1301  ! normal component
1302  elemg = 0.d0
1303  do j=1,nn
1304  elemg(:) = elemg(:)+shapefunc(j)*elemdisp(:,j)
1305  enddo
1306  dg(:) = ddisp(3*slave-2:3*slave) - elemg(:)
1307  dgn = dot_product( contact%states(i)%direction, dg )
1308  contact%states(i)%wkdist = contact%states(i)%wkdist-dgn
1309  contact%states(i)%multiplier(1) = contact%states(i)%multiplier(1)- mu*contact%states(i)%wkdist
1310  contact%states(i)%distance = contact%states(i)%distance - dgn
1311  cnt = cnt+1
1312  lgnt(1) = lgnt(1)- contact%states(i)%wkdist
1313 
1314  if( fcoeff==0.d0 ) cycle
1315  ! tangent component
1316  call dispincrematrix( contact%states(i)%lpos(1:2), etype, nn, elemcrd, tangent, &
1317  metric, dispmat )
1318  dxi(1) = dot_product( dispmat(1,1:nn*3+3), edisp(1:nn*3+3) )
1319  dxi(2) = dot_product( dispmat(2,1:nn*3+3), edisp(1:nn*3+3) )
1320  dxy(:) = matmul( metric, dxi )
1321  fric(1:2) = contact%states(i)%multiplier(2:3) + mut*dxy(1:2)
1322  f3(:) = fric(1)*dispmat(1,:)+fric(2)*dispmat(2,:)
1323  dgn = dsqrt( f3(1)*f3(1)+f3(2)*f3(2)+f3(3)*f3(3) )
1324  if( contact%states(i)%multiplier(1)>0.d0 ) then
1325  if( dgn > fcoeff*contact%states(i)%multiplier(1) ) then
1326  if( contact%states(i)%state==contactstick ) then
1327  ctchanged= .true.
1328  print *, "Node", slave, "to slip state",dgn, fcoeff*contact%states(i)%multiplier(1)
1329  endif
1330  contact%states(i)%state = contactslip
1331  fric(:) = fric(:)*fcoeff*contact%states(i)%multiplier(1)/dgn
1332  else
1333  if( contact%states(i)%state==contactslip ) then
1334  ctchanged= .true.
1335  print *, "Node", slave, "to stick state",dgn, fcoeff*contact%states(i)%multiplier(1)
1336  endif
1337  contact%states(i)%state = contactstick
1338  endif
1339  endif
1340  contact%states(i)%multiplier(2:3) = fric(:)
1341  dxy(:) = matmul( dg, tangent )
1342  lgnt(2) = lgnt(2)+dsqrt( dxy(1)*dxy(1)+dxy(2)*dxy(2) )
1343  enddo
1344  if(cnt>0) lgnt(:) = lgnt(:)/cnt
1345  gnt = gnt + lgnt
1346  end subroutine update_contact_multiplier
1347 
1350  subroutine update_tied_multiplier( contact, disp, ddisp, mu, ctchanged )
1351  type( tcontact ), intent(inout) :: contact
1352  real(kind=kreal), intent(in) :: disp(:)
1353  real(kind=kreal), intent(in) :: ddisp(:)
1354  real(kind=kreal), intent(in) :: mu
1355  logical, intent(inout) :: ctchanged
1356 
1357  integer(kind=kint) :: slave, etype, master
1358  integer(kind=kint) :: nn, i, j, iss, cnt
1359  real(kind=kreal) :: dg(3), dgmax
1360  real(kind=kreal) :: shapefunc(l_max_surface_node)
1361  real(kind=kreal) :: edisp(3*l_max_elem_node+3)
1362 
1363  do i= 1, size(contact%slave)
1364  if( contact%states(i)%state==contactfree ) cycle ! not in contact
1365  slave = contact%slave(i)
1366  edisp(1:3) = disp(3*slave-2:3*slave)+ddisp(3*slave-2:3*slave)
1367  master = contact%states(i)%surface
1368 
1369  nn = size( contact%master(master)%nodes )
1370  etype = contact%master(master)%etype
1371  do j=1,nn
1372  iss = contact%master(master)%nodes(j)
1373  edisp(3*j+1:3*j+3) = disp(3*iss-2:3*iss)+ddisp(3*iss-2:3*iss)
1374  enddo
1375  call getshapefunc( etype, contact%states(i)%lpos(1:2), shapefunc )
1376 
1377  ! normal component
1378  dg(1:3) = edisp(1:3)
1379  do j=1,nn
1380  dg(1:3) = dg(1:3)-shapefunc(j)*edisp(3*j+1:3*j+3)
1381  enddo
1382 
1383  contact%states(i)%multiplier(1:3) = contact%states(i)%multiplier(1:3) + mu*dg(1:3)
1384 
1385  ! check if tied constraint converged
1386  dgmax = 0.d0
1387  do j=1,(nn+1)*3
1388  dgmax = dgmax + dabs(edisp(j))
1389  enddo
1390  dgmax = dgmax/dble((nn+1)*3)
1391  do j=1,3
1392  if( dabs(dg(j))/dmax1(1.d0,dgmax) > 1.d-3 ) ctchanged = .true.
1393  enddo
1394 
1395  enddo
1396  end subroutine
1397 
1399  subroutine ass_contact_force( contact, coord, disp, B )
1400  type( tcontact ), intent(in) :: contact
1401  real(kind=kreal), intent(in) :: coord(:)
1402  real(kind=kreal), intent(in) :: disp(:)
1403  real(kind=kreal), intent(inout) :: b(:)
1404 
1405  integer(kind=kint) :: slave, etype, master
1406  integer(kind=kint) :: nn, i, j, iss
1407  real(kind=kreal) :: fcoeff, nrlforce, tangent(3,2)
1408  real(kind=kreal) :: elemcrd(3, l_max_elem_node )
1409  real(kind=kreal) :: shapefunc(l_max_surface_node)
1410  real(kind=kreal) :: metric(2,2), dispmat(2,l_max_elem_node*3+3)
1411  real(kind=kreal) :: fric(2), f3(l_max_elem_node*3+3)
1412  fcoeff = contact%fcoeff
1413  do i= 1, size(contact%slave)
1414  if( contact%states(i)%state==contactfree ) cycle ! not in contact
1415  slave = contact%slave(i)
1416  master = contact%states(i)%surface
1417 
1418  nn = size( contact%master(master)%nodes )
1419  etype = contact%master(master)%etype
1420  call getshapefunc( etype, contact%states(i)%lpos(1:2), shapefunc )
1421 
1422  ! normal component
1423  nrlforce = contact%states(i)%multiplier(1)
1424  b(3*slave-2:3*slave) = b(3*slave-2:3*slave)-nrlforce*contact%states(i)%direction
1425 
1426  do j=1,nn
1427  iss = contact%master(master)%nodes(j)
1428  b(3*iss-2:3*iss) = b(3*iss-2:3*iss)+nrlforce* &
1429  shapefunc(j)*contact%states(i)%direction
1430  enddo
1431 
1432  if( fcoeff==0.d0 ) cycle
1433  ! tangent component
1434  do j=1,nn
1435  iss = contact%master(master)%nodes(j)
1436  elemcrd(:,j) = coord(3*iss-2:3*iss)+disp(3*iss-2:3*iss)
1437  enddo
1438  call dispincrematrix( contact%states(i)%lpos(1:2), etype, nn, elemcrd, tangent, &
1439  metric, dispmat )
1440 
1441  fric(1:2) = contact%states(i)%multiplier(2:3)
1442  f3(:) = fric(1)*dispmat(1,:)+fric(2)*dispmat(2,:)
1443  b(3*slave-2:3*slave) = b(3*slave-2:3*slave)-f3(1:3)
1444  do j=1,nn
1445  iss = contact%master(master)%nodes(j)
1446  b(3*iss-2:3*iss) = b(3*iss-2:3*iss)-f3(3*j+1:3*j+3)
1447  enddo
1448  enddo
1449 
1450  end subroutine ass_contact_force
1451 
1453  subroutine set_contact_state_vector( contact, dt, relvel_vec, state_vec )
1454  type( tcontact ), intent(in) :: contact
1455  real(kind=kreal), intent(in) :: dt
1456  real(kind=kreal), intent(inout) :: relvel_vec(:)
1457  real(kind=kreal), intent(inout) :: state_vec(:)
1458 
1459  integer(kind=kint) :: i, slave
1460 
1461  do i= 1, size(contact%slave)
1462  slave = contact%slave(i)
1463  if( state_vec(slave) < 0.1d0 .or. contact%states(i)%state > 0 ) &
1464  & state_vec(slave) = dble(contact%states(i)%state)
1465 
1466  if( contact%states(i)%state==contactfree ) cycle ! not in contact
1467  if( dt < 1.d-16 ) cycle ! too small delta t
1468  relvel_vec(3*slave-2:3*slave) = contact%states(i)%reldisp(1:3)/dt
1469  enddo
1470 
1471  end subroutine set_contact_state_vector
1472 
1473  subroutine update_contact_tangentforce( contact )
1474  type( tcontact ), intent(inout) :: contact
1475 
1476  integer(kind=kint) :: i
1477 
1478  do i= 1, size(contact%slave)
1479  if( contact%states(i)%state==contactfree ) then
1480  contact%states(i)%tangentForce(1:3) = 0.d0
1481  contact%states(i)%tangentForce_trial(1:3) = 0.d0
1482  contact%states(i)%tangentForce_final(1:3) = 0.d0
1483  else
1484  contact%states(i)%tangentForce(1:3) = contact%states(i)%tangentForce_final(1:3)
1485  end if
1486  contact%states(i)%tangentForce1(1:3) = contact%states(i)%tangentForce(1:3)
1487  enddo
1488  end subroutine update_contact_tangentforce
1489 
1491  subroutine track_contact_position_exp( nslave, contact, currpos, currdisp, infoCTChange, nodeID, elemID )
1492  integer, intent(in) :: nslave
1493  type( tcontact ), intent(inout) :: contact
1494  type( fstr_info_contactchange ), intent(inout) :: infoCTChange
1495  real(kind=kreal), intent(in) :: currpos(:)
1496  real(kind=kreal), intent(in) :: currdisp(:)
1497  integer(kind=kint), intent(in) :: nodeID(:)
1498  integer(kind=kint), intent(in) :: elemID(:)
1499 
1500  integer(kind=kint) :: slave, sid0, sid, etype
1501  integer(kind=kint) :: nn, i, j, iSS
1502  real(kind=kreal) :: coord(3), elem(3, l_max_elem_node ), elem0(3, l_max_elem_node )
1503  logical :: isin
1504  real(kind=kreal) :: opos(2), odirec(3)
1505  integer(kind=kint) :: bktID, nCand, idm
1506  integer(kind=kint), allocatable :: indexCand(:)
1507 
1508  sid = 0
1509 
1510  slave = contact%slave(nslave)
1511  coord(:) = currpos(3*slave-2:3*slave)
1513  sid0 = contact%states(nslave)%surface
1514  opos = contact%states(nslave)%lpos(1:2)
1515  odirec = contact%states(nslave)%direction
1516  etype = contact%master(sid0)%etype
1517  nn = getnumberofnodes( etype )
1518  do j=1,nn
1519  iss = contact%master(sid0)%nodes(j)
1520  elem(1:3,j)=currpos(3*iss-2:3*iss)
1521  elem0(1:3,j)=currpos(3*iss-2:3*iss)-currdisp(3*iss-2:3*iss)
1522  enddo
1523  call project_point2element( coord,etype,nn,elem,contact%master(sid0)%reflen,contact%states(nslave), &
1524  isin,contact%cparam%DISTCLR_NOCHECK,contact%states(nslave)%lpos(1:2),contact%cparam%CLR_SAME_ELEM )
1525  if( .not. isin ) then
1526  do i=1, contact%master(sid0)%n_neighbor
1527  sid = contact%master(sid0)%neighbor(i)
1528  etype = contact%master(sid)%etype
1529  nn = getnumberofnodes( etype )
1530  do j=1,nn
1531  iss = contact%master(sid)%nodes(j)
1532  elem(1:3,j)=currpos(3*iss-2:3*iss)
1533  enddo
1534  call project_point2element( coord,etype,nn,elem,contact%master(sid)%reflen,contact%states(nslave), &
1535  isin,contact%cparam%DISTCLR_NOCHECK,localclr=contact%cparam%CLEARANCE )
1536  if( isin ) then
1537  contact%states(nslave)%surface = sid
1538  exit
1539  endif
1540  enddo
1541  endif
1542 
1543  if( .not. isin ) then ! such case is considered to rarely or never occur
1544  write(*,*) 'Warning: contact moved beyond neighbor elements'
1545  ! get master candidates from bucketDB
1546  bktid = bucketdb_getbucketid(contact%master_bktDB, coord)
1547  ncand = bucketdb_getnumcand(contact%master_bktDB, bktid)
1548  if (ncand > 0) then
1549  allocate(indexcand(ncand))
1550  call bucketdb_getcand(contact%master_bktDB, bktid, ncand, indexcand)
1551  do idm= 1, ncand
1552  sid = indexcand(idm)
1553  if( sid==sid0 ) cycle
1554  if( any(sid==contact%master(sid0)%neighbor(:)) ) cycle
1555  if (.not. is_in_surface_box( contact%master(sid), coord(1:3), contact%cparam%BOX_EXP_RATE )) cycle
1556  etype = contact%master(sid)%etype
1557  nn = size( contact%master(sid)%nodes )
1558  do j=1,nn
1559  iss = contact%master(sid)%nodes(j)
1560  elem(1:3,j)=currpos(3*iss-2:3*iss)
1561  enddo
1562  call project_point2element( coord,etype,nn,elem,contact%master(sid)%reflen,contact%states(nslave), &
1563  isin,contact%cparam%DISTCLR_NOCHECK,localclr=contact%cparam%CLEARANCE )
1564  if( isin ) then
1565  contact%states(nslave)%surface = sid
1566  exit
1567  endif
1568  enddo
1569  deallocate(indexcand)
1570  endif
1571  endif
1572 
1573  if( isin ) then
1574  if( contact%states(nslave)%surface==sid0 ) then
1575  if(any(dabs(contact%states(nslave)%lpos(1:2)-opos(:)) >= contact%cparam%CLR_DIFFLPOS)) then
1576  !$omp atomic
1577  infoctchange%contact2difflpos = infoctchange%contact2difflpos + 1
1578  endif
1579  else
1580  write(*,'(A,i10,A,i10,A,f7.3,A,2f7.3)') "Node",nodeid(slave)," move to contact with", &
1581  elemid(contact%master(sid)%eid), " with distance ", &
1582  contact%states(nslave)%distance," at ",contact%states(nslave)%lpos(1:2)
1583  !$omp atomic
1584  infoctchange%contact2neighbor = infoctchange%contact2neighbor + 1
1585  endif
1586  iss = isinsideelement( etype, contact%states(nslave)%lpos(1:2), contact%cparam%CLR_CAL_NORM )
1587  if( iss>0 ) then
1588  call cal_node_normal( contact%states(nslave)%surface, iss, contact%master, currpos, &
1589  contact%states(nslave)%lpos(1:2), contact%states(nslave)%direction(:) )
1590  endif
1591  else if( .not. isin ) then
1592  write(*,'(A,i10,A)') "Node",nodeid(slave)," move out of contact"
1593  contact%states(nslave)%state = contactfree
1594  contact%states(nslave)%multiplier(:) = 0.d0
1595  endif
1596 
1597  end subroutine track_contact_position_exp
1598 
1602  subroutine scan_contact_state_exp( contact, currpos, currdisp, infoCTChange, &
1603  nodeID, elemID, is_init, active )
1604  type( tcontact ), intent(inout) :: contact
1605  type( fstr_info_contactchange ), intent(inout) :: infoCTChange
1606  real(kind=kreal), intent(in) :: currpos(:)
1607  real(kind=kreal), intent(in) :: currdisp(:)
1608  integer(kind=kint), intent(in) :: nodeID(:)
1609  integer(kind=kint), intent(in) :: elemID(:)
1610  logical, intent(in) :: is_init
1611  logical, intent(out) :: active
1612 
1613  real(kind=kreal) :: distclr
1614  integer(kind=kint) :: slave, id, etype
1615  integer(kind=kint) :: nn, i, j, iSS, nactive
1616  real(kind=kreal) :: coord(3), elem(3, l_max_elem_node )
1617  real(kind=kreal) :: nlforce
1618  logical :: isin
1619  integer(kind=kint), allocatable :: contact_surf(:), states_prev(:)
1620  !
1621  integer, pointer :: indexMaster(:),indexCand(:)
1622  integer :: nMaster,idm,nMasterMax,bktID,nCand
1623  logical :: is_cand
1624 
1625  if( is_init ) then
1626  distclr = contact%cparam%DISTCLR_INIT
1627  else
1628  distclr = contact%cparam%DISTCLR_FREE
1629  if( contact%algtype == contacttied ) then
1630  active = .false.
1631  do i= 1, size(contact%slave)
1632  if( contact%states(i)%state==contactstick ) then
1633  active = .true.
1634  exit
1635  endif
1636  enddo
1637  return
1638  endif
1639  endif
1640 
1641  allocate(contact_surf(size(nodeid)))
1642  allocate(states_prev(size(contact%slave)))
1643  contact_surf(:) = size(elemid)+1
1644  do i = 1, size(contact%slave)
1645  states_prev(i) = contact%states(i)%state
1646  enddo
1647 
1648  call update_surface_box_info( contact%master, currpos )
1649  call update_surface_bucket_info( contact%master, contact%master_bktDB )
1650 
1651  !$omp parallel do &
1652  !$omp& default(none) &
1653  !$omp& private(i,slave,id,nlforce,coord,indexMaster,nMaster,nn,j,iSS,elem,is_cand,idm,etype,isin, &
1654  !$omp& bktID,nCand,indexCand) &
1655  !$omp& firstprivate(nMasterMax) &
1656  !$omp& shared(contact,infoCTChange,currpos,currdisp,nodeID,elemID,distclr,contact_surf) &
1657  !$omp& reduction(.or.:active) &
1658  !$omp& schedule(dynamic,1)
1659  do i= 1, size(contact%slave)
1660  slave = contact%slave(i)
1661  if( contact%states(i)%state==contactstick .or. contact%states(i)%state==contactslip ) then
1662  call track_contact_position_exp( i, contact, currpos, currdisp, infoctchange, nodeid, elemid )
1663  if( contact%states(i)%state /= contactfree ) then
1664  contact_surf(contact%slave(i)) = -contact%states(i)%surface
1665  endif
1666  else if( contact%states(i)%state==contactfree ) then
1667  coord(:) = currpos(3*slave-2:3*slave)
1668 
1669  ! get master candidates from bucketDB
1670  bktid = bucketdb_getbucketid(contact%master_bktDB, coord)
1671  ncand = bucketdb_getnumcand(contact%master_bktDB, bktid)
1672  if (ncand == 0) cycle
1673  allocate(indexcand(ncand))
1674  call bucketdb_getcand(contact%master_bktDB, bktid, ncand, indexcand)
1675 
1676  nmastermax = ncand
1677  allocate(indexmaster(nmastermax))
1678  nmaster = 0
1679 
1680  ! narrow down candidates
1681  do idm= 1, ncand
1682  id = indexcand(idm)
1683  if (.not. is_in_surface_box( contact%master(id), coord(1:3), contact%cparam%BOX_EXP_RATE )) cycle
1684  nmaster = nmaster + 1
1685  indexmaster(nmaster) = id
1686  enddo
1687  deallocate(indexcand)
1688 
1689  if(nmaster == 0) then
1690  deallocate(indexmaster)
1691  cycle
1692  endif
1693 
1694  do idm = 1,nmaster
1695  id = indexmaster(idm)
1696  etype = contact%master(id)%etype
1697  nn = size( contact%master(id)%nodes )
1698  do j=1,nn
1699  iss = contact%master(id)%nodes(j)
1700  elem(1:3,j)=currpos(3*iss-2:3*iss)
1701  enddo
1702  call project_point2element( coord,etype,nn,elem,contact%master(id)%reflen,contact%states(i), &
1703  isin,distclr,localclr=contact%cparam%CLEARANCE )
1704  if( .not. isin ) cycle
1705  contact%states(i)%surface = id
1706  contact%states(i)%multiplier(:) = 0.d0
1707  iss = isinsideelement( etype, contact%states(i)%lpos(1:2), contact%cparam%CLR_CAL_NORM )
1708  if( iss>0 ) &
1709  call cal_node_normal( id, iss, contact%master, currpos, contact%states(i)%lpos(1:2), &
1710  contact%states(i)%direction(:) )
1711  contact_surf(contact%slave(i)) = id
1712  write(*,'(A,i10,A,i10,A,f7.3,A,2f7.3,A,3f7.3)') "Node",nodeid(slave)," contact with element", &
1713  elemid(contact%master(id)%eid), &
1714  " with distance ", contact%states(i)%distance," at ",contact%states(i)%lpos(1:2), &
1715  " along direction ", contact%states(i)%direction
1716  exit
1717  enddo
1718  deallocate(indexmaster)
1719  endif
1720  enddo
1721  !$omp end parallel do
1722 
1723  call hecmw_contact_comm_allreduce_i(contact%comm, contact_surf, hecmw_min)
1724  nactive = 0
1725  do i = 1, size(contact%slave)
1726  if (contact%states(i)%state /= contactfree) then ! any slave in contact
1727  if (abs(contact_surf(contact%slave(i))) /= contact%states(i)%surface) then ! that is in contact with other surface
1728  contact%states(i)%state = contactfree ! should be freed
1729  write(*,'(A,i10,A,i6,A,i6,A)') "Node",nodeid(contact%slave(i)), &
1730  " in rank",hecmw_comm_get_rank()," freed due to duplication"
1731  else
1732  nactive = nactive + 1
1733  endif
1734  endif
1735  if (states_prev(i) == contactfree .and. contact%states(i)%state /= contactfree) then
1736  infoctchange%free2contact = infoctchange%free2contact + 1
1737  elseif (states_prev(i) /= contactfree .and. contact%states(i)%state == contactfree) then
1738  infoctchange%contact2free = infoctchange%contact2free + 1
1739  endif
1740  enddo
1741  active = (nactive > 0)
1742  deallocate(contact_surf)
1743  deallocate(states_prev)
1744  end subroutine scan_contact_state_exp
1745 
1746 end module mcontactdef
This module provides bucket-search functionality It provides definition of bucket info and its access...
integer(kind=kint) function, public bucketdb_getbucketid(bktdb, x)
Get bucket ID that includes given point.
subroutine, public bucketdb_finalize(bktdb)
Finalizer.
integer(kind=kint) function, public bucketdb_getnumcand(bktdb, bid)
Get number of candidates within neighboring buckets of a given bucket Bucket ID has to be obtained wi...
subroutine, public bucketdb_init(bktdb)
Initializer.
subroutine, public bucketdb_getcand(bktdb, bid, ncand, cand)
Get candidates within neighboring buckets of a given bucket Number of candidates has to be obtained w...
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
subroutine getvertexcoord(fetype, cnode, localcoord)
Get the natural coord of a vertex node.
Definition: element.f90:1143
real(kind=kreal) function, dimension(3) surfacenormal(fetype, nn, localcoord, elecoord)
Calculate normal of 3d-surface.
Definition: element.f90:870
Definition: hecmw.f90:6
This module provide functions of contact stiffness calculation.
Definition: contact_lib.f90:6
integer, parameter contacttied
contact type or algorithm definition
Definition: contact_lib.f90:21
integer, parameter contactfslid
Definition: contact_lib.f90:24
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.
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.
subroutine update_tangentforce(etype, nn, elemt0, elemt, cstate)
This subroutine find the projection of a slave point onto master surface.
subroutine contact_state_init(cstate)
Initializer.
Definition: contact_lib.f90:57
subroutine set_shrink_factor(ctime, cstate, etime, if_type)
integer, parameter contactslip
Definition: contact_lib.f90:18
integer, parameter contactfree
contact state definition
Definition: contact_lib.f90:16
integer, parameter c_if_slave
contact interference type
Definition: contact_lib.f90:27
subroutine dispincrematrix(pos, etype, nnode, ele, tangent, tensor, matrix)
This subroutine calculate the relation between global disp and displacement along natural coordinate ...
integer, parameter contactstick
Definition: contact_lib.f90:17
subroutine, public hecmw_contact_comm_allreduce_i(conComm, vec, op)
subroutine, public hecmw_contact_comm_init(conComm, hecMESH, ndof, n_contact_dof, contact_dofs)
subroutine, public hecmw_contact_comm_finalize(conComm)
This module manages the data structure for contact calculation.
subroutine set_contact_state_vector(contact, dt, relvel_vec, state_vec)
This subroutine setup contact output nodal vectors.
subroutine update_contact_multiplier(contact, coord, disp, ddisp, fcoeff, mu, mut, gnt, ctchanged)
This subroutine update lagrangian multiplier and the distance between contacting nodes.
integer(kind=kint) function check_apply_contact_if(contact_if, contacts)
subroutine update_tied_multiplier(contact, disp, ddisp, mu, ctchanged)
This subroutine update lagrangian multiplier and the distance between contacting nodes.
logical function fstr_contact_check(contact, hecMESH)
Check the consistency with given mesh of contact definition.
subroutine clear_contact_state(contact)
Reset contact state all to free.
subroutine update_contact_tangentforce(contact)
subroutine track_contact_position_exp(nslave, contact, currpos, currdisp, infoCTChange, nodeID, elemID)
This subroutine tracks down next contact position after a finite slide.
subroutine scan_contact_state_exp(contact, currpos, currdisp, infoCTChange, nodeID, elemID, is_init, active)
This subroutine update contact states, which include.
subroutine scan_embed_state(flag_ctAlgo, embed, currpos, currdisp, ndforce, infoCTChange, nodeID, elemID, is_init, active, mu, B)
This subroutine update contact states, which include.
subroutine fstr_write_contact(file, contact)
Write out contact definition.
subroutine fstr_contact_finalize(contact)
Finalizer.
subroutine update_direction(nslave, contact, currpos)
This subroutine tracks down next contact position after a finite slide.
logical function fstr_embed_init(embed, hecMESH, cparam, myrank)
Initializer of tContactState for embed case.
subroutine ass_contact_force(contact, coord, disp, B)
This subroutine assemble contact force into contacing nodes.
subroutine calcu_contact_force0(contact, coord, disp, ddisp, fcoeff, mu, mut, B)
This subroutine updates contact condition as follows:
subroutine scan_contact_state(flag_ctAlgo, contact, currpos, currdisp, ndforce, infoCTChange, nodeID, elemID, is_init, active, mu, B)
This subroutine update contact states, which include.
logical function fstr_contact_init(contact, hecMESH, cparam, myrank)
Initializer of tContactState.
subroutine calcu_tied_force0(contact, disp, ddisp, mu, B)
This subroutine updates contact condition as follows:
This module manage the parameters for contact calculation.
This module manages surface elements in 3D It provides basic definition of surface elements (triangla...
Definition: surf_ele.f90:8
subroutine initialize_surf(eid, etype, nsurf, surf)
Initializer.
Definition: surf_ele.f90:38
subroutine update_surface_reflen(surf, coord)
Compute reference length of surface elements.
Definition: surf_ele.f90:220
subroutine write_surf(file, surf)
Write out elemental surface.
Definition: surf_ele.f90:76
subroutine update_surface_box_info(surf, currpos)
Update info of cubic box including surface elements.
Definition: surf_ele.f90:240
subroutine find_surface_neighbor(surf, bktDB)
Find neighboring surface elements.
Definition: surf_ele.f90:88
subroutine update_surface_bucket_info(surf, bktDB)
Update bucket info for searching surface elements.
Definition: surf_ele.f90:280
subroutine finalize_surf(surf)
Memory management subroutine.
Definition: surf_ele.f90:69
logical function is_in_surface_box(surf, x0, exp_rate)
Check if given point is within cubic box including surface element.
Definition: surf_ele.f90:265
This structure records contact status.
Definition: contact_lib.f90:30
Structure to includes all info needed by contact calculation.
Structure to define surface group.
Definition: surf_ele.f90:19