FrontISTR  5.9.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 elementinfo
16  use msurfelement
18  use bucket_search
19  use mcontactparam
20 
21  implicit none
22 
23  real(kind=kreal), save :: cgn=1.d-5
24  real(kind=kreal), save :: cgt=1.d-3
25 
26  real(kind=kreal), save :: gnt(2)
28  real(kind=kreal), save :: bakgnt(2)
30 
31  integer, parameter :: contactunknown = -1
33  integer, parameter :: contactfree = -1
34  integer, parameter :: contactstick = 1
35  integer, parameter :: contactslip = 2
36 
38  integer, parameter :: contacttied = 1
39  integer, parameter :: contactglued = 2
40  integer, parameter :: contactsslid = 3
41  integer, parameter :: contactfslid = 4
42 
44  integer, parameter :: c_if_slave = 1
45  integer, parameter :: c_if_master = 2
46 
48  integer, parameter :: kctforresidual = 1
49  integer, parameter :: kctforoutput = 2
50 
53  integer :: state
54  integer :: surface
55  real(kind=kreal) :: distance
56  real(kind=kreal) :: wkdist
57  real(kind=kreal) :: lpos(3)
58  real(kind=kreal) :: gpos(3)
59  real(kind=kreal) :: direction(3)
60  real(kind=kreal) :: multiplier(3)
62  real(kind=kreal) :: tangentforce(3)
63  real(kind=kreal) :: tangentforce1(3)
64  real(kind=kreal) :: tangentforce_trial(3)
65  real(kind=kreal) :: tangentforce_final(3)
66  real(kind=kreal) :: reldisp(3)
67  !
68  real(kind=kreal) :: shrink_factor
69  real(kind=kreal) :: time_factor
70  real(kind=kreal) :: init_pos
71  real(kind=kreal) :: end_pos
72  integer :: interference_flag
73  end type
74 
76  type tcontact
77  ! following contact definition
78  character(len=HECMW_NAME_LEN) :: name
79  integer :: ctype
80  integer :: group
81  character(len=HECMW_NAME_LEN) :: pair_name
82  integer :: surf_id1, surf_id2
83  integer :: surf_id1_sgrp
84  type(tsurfelement), pointer :: master(:)=>null()
85  integer, pointer :: slave(:)=>null()
86  real(kind=kreal) :: fcoeff
87  real(kind=kreal) :: npenalty
88  real(kind=kreal) :: tpenalty
89  real(kind=kreal) :: refstiff
90 
91  real(kind=kreal) :: ctime
92  integer(kind=kint) :: if_type
93  real(kind=kreal) :: if_etime
94  real(kind=kreal) :: initial_pos
95  real(kind=kreal) :: end_pos
96  ! following algorithm
97  ! -1: not initialized
98  ! 1: TIED-Just rigidly fixed the two surfaces
99  ! 2: GLUED-Distance between the two surfaces to zero and glue them
100  ! 3: SSLID-Small sliding contact( no position but with contact state change)
101  ! 4: FSLID-Finite sliding contact (both changes in contact state and position possible)
102  integer :: algtype
103 
104  logical :: mpced
105  logical :: symmetric
106 
107  ! following contact state
108  type(tcontactstate), pointer :: states(:)=>null()
109 
110  type(hecmwst_contact_comm) :: comm
111  type(bucketdb) :: master_bktdb
112 
113  type(tcontactparam), pointer :: cparam=>null()
114  end type tcontact
115 
117  logical :: active
118  integer(kind=kint) :: contact2free
119  integer(kind=kint) :: contact2neighbor
120  integer(kind=kint) :: contact2difflpos
121  integer(kind=kint) :: free2contact
122  integer(kind=kint) :: contactnode_previous
123  integer(kind=kint) :: contactnode_current
124  end type fstr_info_contactchange
125 
126  private :: is_mpc_available
127  private :: is_active_contact
128 
129 contains
130 
132  subroutine contact_state_init(cstate)
133  type(tcontactstate), intent(inout) :: cstate
134  cstate%state = -1
135  cstate%surface = -1
136  cstate%distance = 0.0d0
137  cstate%wkdist = 0.0d0
138  cstate%lpos(:) = 0.0d0
139  cstate%gpos(:) = 0.0d0
140  cstate%direction(:) = 0.0d0
141  cstate%multiplier(:) = 0.0d0
142  cstate%tangentForce(:) = 0.0d0
143  cstate%tangentForce1(:) = 0.0d0
144  cstate%tangentForce_trial(:) = 0.0d0
145  cstate%tangentForce_final(:) = 0.0d0
146  cstate%reldisp(:) = 0.0d0
147  cstate%shrink_factor = 0.0d0
148  cstate%time_factor = 0.0d0
149  cstate%init_pos = 0.0d0
150  cstate%end_pos = 0.0d0
151  cstate%interference_flag = 0
152  end subroutine
153 
155  subroutine contact_state_copy(cstate1, cstate2)
156  type(tcontactstate), intent(in) :: cstate1
157  type(tcontactstate), intent(inout) :: cstate2
158  cstate2 = cstate1
159  end subroutine
160 
162  subroutine print_contact_state(fnum, cstate)
163  integer, intent(in) :: fnum
164  type(tcontactstate), intent(in) :: cstate
165  write(fnum, *) "--Contact state=",cstate%state
166  write(fnum, *) cstate%surface, cstate%distance
167  write(fnum, *) cstate%lpos
168  write(fnum, *) cstate%direction
169  write(fnum, *) cstate%multiplier
170  end subroutine
171 
173  logical function is_mpc_available( contact )
174  type(tcontact), intent(in) :: contact
175  is_mpc_available = .true.
176  if( contact%fcoeff/=0.d0 ) is_mpc_available = .false.
177  end function
178 
180  subroutine fstr_write_contact( file, contact )
181  integer(kind=kint), intent(in) :: file
182  type(tcontact), intent(in) :: contact
183  integer :: i
184  write(file,*) "CONTACT:", contact%ctype,contact%group,trim(contact%pair_name),contact%fcoeff
185  write(file,*) "---Slave----"
186  write(file,*) "num.slave",size(contact%slave)
187  if( associated(contact%slave) ) then
188  do i=1,size(contact%slave)
189  write(file, *) contact%slave(i)
190  enddo
191  endif
192  write(file,*) "----master---"
193  write(file,*) "num.master",size(contact%master)
194  if( associated(contact%master) ) then
195  do i=1,size(contact%master)
196  call write_surf( file, contact%master(i) )
197  enddo
198  endif
199  end subroutine
200 
202  subroutine fstr_contact_finalize( contact )
203  type(tcontact), intent(inout) :: contact
204  integer :: i
205  if( associated( contact%slave ) ) deallocate(contact%slave)
206  if( associated( contact%master ) ) then
207  do i=1,size( contact%master )
208  call finalize_surf( contact%master(i) )
209  enddo
210  deallocate(contact%master)
211  endif
212  if( associated(contact%states) ) deallocate(contact%states)
213  call hecmw_contact_comm_finalize(contact%comm)
214  call bucketdb_finalize( contact%master_bktDB )
215  end subroutine
216 
218  logical function fstr_contact_check( contact, hecMESH )
219  type(tcontact), intent(inout) :: contact
220  type(hecmwst_local_mesh), pointer :: hecmesh
221 
222  integer :: i
223  logical :: isfind
224 
225  fstr_contact_check = .false.
226 
227  ! if contact pair exist?
228  isfind = .false.
229  do i=1,hecmesh%contact_pair%n_pair
230  if( hecmesh%contact_pair%name(i) == contact%pair_name ) then
231  contact%ctype = hecmesh%contact_pair%type(i)
232  contact%surf_id1 = hecmesh%contact_pair%slave_grp_id(i)
233  contact%surf_id2 = hecmesh%contact_pair%master_grp_id(i)
234  contact%surf_id1_sgrp = hecmesh%contact_pair%slave_orisgrp_id(i)
235  isfind = .true.
236  endif
237  enddo
238  if( .not. isfind ) return;
239  if( contact%fcoeff<=0.d0 ) contact%fcoeff=0.d0
240  if( contact%ctype < 1 .and. contact%ctype > 3 ) return
241  if( contact%group<=0 ) return
242 
243  fstr_contact_check = .true.
244  end function
245 
247  logical function fstr_contact_init( contact, hecMESH, cparam, myrank )
248  type(tcontact), intent(inout) :: contact
249  type(hecmwst_local_mesh), pointer :: hecmesh
250  type(tcontactparam), target :: cparam
251  integer(kint),intent(in),optional :: myrank
252 
253  integer :: i, j, is, ie, cgrp, nsurf, nslave, ic, ic_type, iss, nn, ii
254  integer :: count,nodeid
255 
256  fstr_contact_init = .false.
257 
258  contact%cparam => cparam
259 
260  ! master surface
261  cgrp = contact%surf_id2
262  if( cgrp<=0 ) return
263  is= hecmesh%surf_group%grp_index(cgrp-1) + 1
264  ie= hecmesh%surf_group%grp_index(cgrp )
265 
266  if(present(myrank)) then
267  ! PARA_CONTACT
268  count = 0
269  do i=is,ie
270  ic = hecmesh%surf_group%grp_item(2*i-1)
271  if(hecmesh%elem_ID(ic*2) /= myrank) cycle
272  count = count + 1
273  enddo
274  allocate( contact%master(count) )
275  count = 0
276  do i=is,ie
277  ic = hecmesh%surf_group%grp_item(2*i-1)
278  if(hecmesh%elem_ID(ic*2) /= myrank) cycle
279  count = count + 1
280  nsurf = hecmesh%surf_group%grp_item(2*i)
281  ic_type = hecmesh%elem_type(ic)
282  call initialize_surf( ic, ic_type, nsurf, contact%master(count) )
283  iss = hecmesh%elem_node_index(ic-1)
284  do j=1, size( contact%master(count)%nodes )
285  nn = contact%master(count)%nodes(j)
286  contact%master(count)%nodes(j) = hecmesh%elem_node_item( iss+nn )
287  enddo
288  enddo
289 
290  else
291  ! not PARA_CONTACT
292  allocate( contact%master(ie-is+1) )
293  do i=is,ie
294  ic = hecmesh%surf_group%grp_item(2*i-1)
295  nsurf = hecmesh%surf_group%grp_item(2*i)
296  ic_type = hecmesh%elem_type(ic)
297  call initialize_surf( ic, ic_type, nsurf, contact%master(i-is+1) )
298  iss = hecmesh%elem_node_index(ic-1)
299  do j=1, size( contact%master(i-is+1)%nodes )
300  nn = contact%master(i-is+1)%nodes(j)
301  contact%master(i-is+1)%nodes(j) = hecmesh%elem_node_item( iss+nn )
302  enddo
303  enddo
304 
305  endif
306 
307  call update_surface_reflen( contact%master, hecmesh%node )
308 
309  cgrp = contact%surf_id1
310  if( cgrp<=0 ) return
311  is= hecmesh%node_group%grp_index(cgrp-1) + 1
312  ie= hecmesh%node_group%grp_index(cgrp )
313  nslave = 0
314  do i=is,ie
315  nodeid = hecmesh%global_node_ID(hecmesh%node_group%grp_item(i))
316  if(present(myrank)) then
317  ! PARA_CONTACT
318  nslave = nslave + 1
319  else
320  ! not PARA_CONTACT
321  if( hecmesh%node_group%grp_item(i) <= hecmesh%nn_internal) then
322  nslave = nslave + 1
323  endif
324  endif
325  enddo
326  allocate( contact%slave(nslave) )
327  ii = 0
328  do i=is,ie
329  if(.not.present(myrank)) then
330  ! not PARA_CONTACT
331  if( hecmesh%node_group%grp_item(i) > hecmesh%nn_internal) cycle
332  endif
333  ii = ii + 1
334  contact%slave(ii) = hecmesh%node_group%grp_item(i)
335  enddo
336 
337  ! contact state
338  allocate( contact%states(nslave) )
339  do i=1,nslave
340  call contact_state_init( contact%states(i) )
341  enddo
342 
343  ! neighborhood of surface group
344  call update_surface_box_info( contact%master, hecmesh%node )
345  call bucketdb_init( contact%master_bktDB )
346  call update_surface_bucket_info( contact%master, contact%master_bktDB )
347  call find_surface_neighbor( contact%master, contact%master_bktDB )
348 
349  ! initialize contact communication table
350  call hecmw_contact_comm_init( contact%comm, hecmesh, 1, nslave, contact%slave )
351 
352  contact%symmetric = .true.
353  fstr_contact_init = .true.
354  end function
355 
357  logical function fstr_embed_init( embed, hecMESH, cparam, myrank )
358  type(tcontact), intent(inout) :: embed
359  type(hecmwst_local_mesh), pointer :: hecmesh
360  type(tcontactparam), target :: cparam
361  integer(kint),intent(in),optional :: myrank
362 
363  integer :: i, j, is, ie, cgrp, nsurf, nslave, ic, ic_type, iss, nn, ii
364  integer :: count,nodeid
365 
366  fstr_embed_init = .false.
367 
368  embed%cparam => cparam
369 
370  ! master surface
371  cgrp = embed%surf_id2
372  if( cgrp<=0 ) return
373  is= hecmesh%elem_group%grp_index(cgrp-1) + 1
374  ie= hecmesh%elem_group%grp_index(cgrp )
375 
376  if(present(myrank)) then
377  ! PARA_CONTACT
378  count = 0
379  do i=is,ie
380  ic = hecmesh%elem_group%grp_item(i)
381  if(hecmesh%elem_ID(ic*2) /= myrank) cycle
382  count = count + 1
383  enddo
384  allocate( embed%master(count) )
385  count = 0
386  do i=is,ie
387  ic = hecmesh%elem_group%grp_item(i)
388  if(hecmesh%elem_ID(ic*2) /= myrank) cycle
389  count = count + 1
390  ic_type = hecmesh%elem_type(ic)
391  call initialize_surf( ic, ic_type, 0, embed%master(count) )
392  iss = hecmesh%elem_node_index(ic-1)
393  do j=1, size( embed%master(count)%nodes )
394  nn = embed%master(count)%nodes(j)
395  embed%master(count)%nodes(j) = hecmesh%elem_node_item( iss+nn )
396  enddo
397  enddo
398 
399  else
400  ! not PARA_CONTACT
401  allocate( embed%master(ie-is+1) )
402  do i=is,ie
403  ic = hecmesh%elem_group%grp_item(i)
404  ic_type = hecmesh%elem_type(ic)
405  call initialize_surf( ic, ic_type, 0, embed%master(i-is+1) )
406  iss = hecmesh%elem_node_index(ic-1)
407  do j=1, size( embed%master(i-is+1)%nodes )
408  nn = embed%master(i-is+1)%nodes(j)
409  embed%master(i-is+1)%nodes(j) = hecmesh%elem_node_item( iss+nn )
410  enddo
411  enddo
412 
413  endif
414 
415  ! slave surface
416  cgrp = embed%surf_id1
417  if( cgrp<=0 ) return
418  is= hecmesh%node_group%grp_index(cgrp-1) + 1
419  ie= hecmesh%node_group%grp_index(cgrp )
420  nslave = 0
421  do i=is,ie
422  nodeid = hecmesh%global_node_ID(hecmesh%node_group%grp_item(i))
423  if(present(myrank)) then
424  ! PARA_CONTACT
425  nslave = nslave + 1
426  else
427  ! not PARA_CONTACT
428  if( hecmesh%node_group%grp_item(i) <= hecmesh%nn_internal) then
429  nslave = nslave + 1
430  endif
431  endif
432  enddo
433  allocate( embed%slave(nslave) )
434  ii = 0
435  do i=is,ie
436  if(.not.present(myrank)) then
437  ! not PARA_CONTACT
438  if( hecmesh%node_group%grp_item(i) > hecmesh%nn_internal) cycle
439  endif
440  ii = ii + 1
441  embed%slave(ii) = hecmesh%node_group%grp_item(i)
442  enddo
443 
444  ! embed state
445  allocate( embed%states(nslave) )
446  do i=1,nslave
447  call contact_state_init( embed%states(i) )
448  enddo
449 
450  ! neighborhood of surface group
451  call update_surface_box_info( embed%master, hecmesh%node )
452  call bucketdb_init( embed%master_bktDB )
453  call update_surface_bucket_info( embed%master, embed%master_bktDB )
454  call find_surface_neighbor( embed%master, embed%master_bktDB )
455 
456  ! initialize contact communication table
457  call hecmw_contact_comm_init( embed%comm, hecmesh, 1, nslave, embed%slave )
458 
459  ! initialize penalty coefficients
460  embed%nPenalty = 1.0d0 ! default normal penalty coefficient
461  embed%tPenalty = 0.1d0 ! default tangential penalty coefficient
462  embed%refStiff = 0.0d0 ! will be calculated after first stiffness assembly
463 
464  embed%symmetric = .true.
465  fstr_embed_init = .true.
466  end function
467 
468  function check_apply_contact_if( contact_if, contacts )
469  type(tcontactinterference), intent(inout) :: contact_if
470  type(tcontact) :: contacts(:)
471 
472  integer :: i, j
473  logical :: isfind
474  integer(kind=kint) :: check_apply_contact_if
475 
477  ! if contact pair exist?
478  isfind = .false.
479  do i = 1, size(contacts)
480  if( contacts(i)%pair_name == contact_if%cp_name ) then
481  contacts(i)%if_type = contact_if%if_type
482  contacts(i)%if_etime = contact_if%etime
483  contacts(i)%initial_pos = contact_if%initial_pos
484  contacts(i)%end_pos = contact_if%end_pos
485  do j = 1, size(contacts(i)%states)
486  contacts(i)%states(j)%interference_flag = contact_if%if_type
487  contacts(i)%states(j)%init_pos = contact_if%initial_pos
488  contacts(i)%states(j)%end_pos = contact_if%end_pos
489  if( contact_if%if_type /= c_if_slave )then
490  contacts(i)%states(j)%time_factor = (contact_if%end_pos - contact_if%initial_pos) / contact_if%etime
491  else
492  contacts(i)%states(j)%time_factor = contact_if%etime
493  end if
494  end do
495  isfind = .true.
496  check_apply_contact_if = 0; return
497  endif
498  enddo
499  if( .not. isfind ) return;
501 
502  end function
503 
505  subroutine clear_contact_state( contact )
506  type(tcontact), intent(inout) :: contact
507  integer :: i
508  if( .not. associated(contact%states) ) return
509  do i=1,size( contact%states )
510  contact%states(i)%state = -1
511  enddo
512  end subroutine
513 
515  logical function is_active_contact( acgrp, contact )
516  integer, intent(in) :: acgrp(:)
517  type(tcontact), intent(in) :: contact
518  if( any( acgrp==contact%group ) ) then
519  is_active_contact = .true.
520  else
521  is_active_contact = .false.
522  endif
523  end function
524 
526  subroutine print_contatct_pair( file, pair )
527  integer(kind=kint), intent(in) :: file
528  type( hecmwst_contact_pair ), intent(in) :: pair
529 
530  integer(kind=kint) :: i
531  write(file,*) "Number of contact pair", pair%n_pair
532  do i=1,pair%n_pair
533  write(file,*) trim(pair%name(i)), pair%type(i), pair%slave_grp_id(i) &
534  ,pair%master_grp_id(i), pair%slave_orisgrp_id(i)
535  enddo
536  end subroutine
537 
538 
539 end module mcontactdef
This module provides bucket-search functionality It provides definition of bucket info and its access...
subroutine, public bucketdb_finalize(bktdb)
Finalizer.
subroutine, public bucketdb_init(bktdb)
Initializer.
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
Definition: hecmw.f90:6
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.
integer, parameter contactslip
subroutine contact_state_copy(cstate1, cstate2)
Copy.
integer, parameter contactsslid
integer, parameter kctforresidual
purpose flag for contact force calculation
integer(kind=kint) function check_apply_contact_if(contact_if, contacts)
integer, parameter contactunknown
real(kind=kreal), dimension(2), save gnt
1:current average penetration; 2:current relative tangent displacement
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.
integer, parameter contactglued
subroutine contact_state_init(cstate)
Initializer.
integer, parameter contacttied
contact type or algorithm definition
integer, parameter kctforoutput
compute contact force for output (CONT_NFORCE/CONT_FRIC)
subroutine fstr_write_contact(file, contact)
Write out contact definition.
subroutine print_contatct_pair(file, pair)
Write out the contact definition read from mesh file.
subroutine print_contact_state(fnum, cstate)
Print out contact state.
subroutine fstr_contact_finalize(contact)
Finalizer.
logical function fstr_embed_init(embed, hecMESH, cparam, myrank)
Initializer of tContactState for embed case.
integer, parameter c_if_slave
contact interference type
real(kind=kreal), dimension(2), save bakgnt
1:current average penetration; 2:current relative tangent displacement!
real(kind=kreal), save cgn
convergent condition of penetration
integer, parameter c_if_master
integer, parameter contactfree
contact state definition
integer, parameter contactfslid
logical function fstr_contact_init(contact, hecMESH, cparam, myrank)
Initializer of tContactState.
real(kind=kreal), save cgt
convergent condition of relative tangent disp
integer, parameter contactstick
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:42
subroutine update_surface_reflen(surf, coord)
Compute reference length of surface elements.
Definition: surf_ele.f90:219
subroutine write_surf(file, surf)
Write out elemental surface.
Definition: surf_ele.f90:79
subroutine update_surface_box_info(surf, currpos)
Update info of cubic box including surface elements.
Definition: surf_ele.f90:238
subroutine find_surface_neighbor(surf, bktDB)
Find neighboring surface elements.
Definition: surf_ele.f90:91
subroutine update_surface_bucket_info(surf, bktDB)
Update bucket info for searching surface elements.
Definition: surf_ele.f90:278
subroutine finalize_surf(surf)
Memory management subroutine.
Definition: surf_ele.f90:72
Structure to includes all info needed by contact calculation.
This structure records contact status.
Structure to define surface group.
Definition: surf_ele.f90:23