FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_contact.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 !-------------------------------------------------------------------------------
33 module mcontact
34 
35  use mcontactdef
36  use hecmw
37  use m_fstr
44  implicit none
45 
46  logical, private :: active
47 
48  ! Re-export functions from sub-modules
50  public :: initialize_embed_vectors
52  public :: calc_contact_area
56 
57 contains
58 
59  subroutine fstr_addcontactstiffness(cstep,ctAlgo, iter,hecMESH,conMAT,hecLagMAT,fstrSOLID)
60 
61  integer(kind=kint) :: cstep
62  integer(kind=kint) :: ctAlgo
63  integer(kind=kint) :: iter
64  type(hecmwst_local_mesh) :: hecMESH
65  type(hecmwst_matrix) :: conMAT
66  type(fstr_solid) :: fstrSOLID
67  type(hecmwst_matrix_lagrange) :: hecLagMAT
68  integer(kind=kint) :: i, grpid
69 
70  if( associated(heclagmat%AL_lagrange) ) heclagmat%AL_lagrange = 0.0d0
71  if( associated(heclagmat%AU_lagrange) ) heclagmat%AU_lagrange = 0.0d0
72 
73  do i = 1, fstrsolid%n_contacts
74 
75  grpid = fstrsolid%contacts(i)%group
76  if( .not. fstr_iscontactactive( fstrsolid, grpid, cstep ) ) cycle
77 
78  call calcu_contact_stiffness_nodesurf( ctalgo, fstrsolid%contacts(i), hecmesh%node(:), fstrsolid%unode(:), &
79  fstrsolid%dunode(:), iter, heclagmat%Lagrange(:), conmat, heclagmat)
80 
81  enddo
82 
83  do i = 1, fstrsolid%n_embeds
84 
85  grpid = fstrsolid%embeds(i)%group
86  if( .not. fstr_isembedactive( fstrsolid, grpid, cstep ) ) cycle
87 
88  call calcu_contact_stiffness_nodesurf( ctalgo, fstrsolid%embeds(i), hecmesh%node(:), fstrsolid%unode(:), &
89  fstrsolid%dunode(:), iter, heclagmat%Lagrange(:), conmat, heclagmat)
90 
91  enddo
92 
93  end subroutine fstr_addcontactstiffness
94 
96  subroutine fstr_update_ndforce_contact(cstep,ctAlgo,hecMESH,hecLagMAT,fstrSOLID,conMAT)
97  integer(kind=kint), intent(in) :: cstep
98  integer(kind=kint), intent(in) :: ctAlgo
99  type(hecmwst_local_mesh) :: hecMESH
100  type(fstr_solid) :: fstrSOLID
101  type(hecmwst_matrix_lagrange) :: hecLagMAT
102  type(hecmwst_matrix) :: conMAT
103 
104  call fstr_contact_ndforce_core(kctforresidual,cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
105 
106  end subroutine fstr_update_ndforce_contact
107 
109  subroutine fstr_calc_contact_output_force(cstep,ctAlgo,hecMESH,hecLagMAT,fstrSOLID,conMAT)
110  integer(kind=kint), intent(in) :: cstep
111  integer(kind=kint), intent(in) :: ctAlgo
112  type(hecmwst_local_mesh) :: hecMESH
113  type(fstr_solid) :: fstrSOLID
114  type(hecmwst_matrix_lagrange) :: hecLagMAT
115  type(hecmwst_matrix) :: conMAT
116 
117  call fstr_contact_ndforce_core(kctforoutput,cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
118 
119  end subroutine fstr_calc_contact_output_force
120 
124  subroutine fstr_contact_ndforce_core(purpose,cstep,ctAlgo,hecMESH,hecLagMAT,fstrSOLID,conMAT)
125  integer(kind=kint), intent(in) :: purpose
126  integer(kind=kint), intent(in) :: cstep
127  integer(kind=kint), intent(in) :: ctAlgo
128  type(hecmwst_local_mesh) :: hecMESH
129  type(fstr_solid) :: fstrSOLID
130  type(hecmwst_matrix_lagrange) :: hecLagMAT
131  type(hecmwst_matrix) :: conMAT
132  integer(kind=kint) :: i, grpid
133 
134  ! Zero-clear output arrays only when computing for output
135  if( purpose == kctforoutput ) then
136  if( associated(fstrsolid%CONT_NFORCE) ) fstrsolid%CONT_NFORCE(:) = 0.d0
137  if( associated(fstrsolid%CONT_FRIC) ) fstrsolid%CONT_FRIC(:) = 0.d0
138  if( associated(fstrsolid%EMBED_NFORCE) ) fstrsolid%EMBED_NFORCE(:) = 0.d0
139  endif
140 
141  do i = 1, fstrsolid%n_contacts
142  grpid = fstrsolid%contacts(i)%group
143  if( .not. fstr_iscontactactive( fstrsolid, grpid, cstep ) ) cycle
144  call calcu_contact_ndforce_nodesurf( purpose, ctalgo, fstrsolid%contacts(i), hecmesh%node(:), fstrsolid%unode(:), &
145  fstrsolid%dunode(:), heclagmat%Lagrange(:), conmat, &
146  fstrsolid%CONT_NFORCE, fstrsolid%CONT_FRIC, heclagmat )
147  enddo
148 
149  do i = 1, fstrsolid%n_embeds
150  grpid = fstrsolid%embeds(i)%group
151  if( .not. fstr_isembedactive( fstrsolid, grpid, cstep ) ) cycle
152  call calcu_contact_ndforce_nodesurf( purpose, ctalgo, fstrsolid%embeds(i), hecmesh%node(:), fstrsolid%unode(:), &
153  fstrsolid%dunode(:), heclagmat%Lagrange(:), conmat, &
154  fstrsolid%EMBED_NFORCE, fstrsolid%EMBED_NFORCE, heclagmat )
155  enddo
156 
157  end subroutine fstr_contact_ndforce_core
158 
160  subroutine fstr_calc_contact_refstiff(cstep, hecMESH, hecMAT, fstrSOLID)
162  integer(kind=kint), intent(in) :: cstep
163  type(hecmwst_local_mesh) :: hecMESH
164  type(hecmwst_matrix) :: hecMAT
165  type(fstr_solid) :: fstrSOLID
166 
167  integer(kind=kint) :: i, grpid, ndof
168  real(kind=kreal), pointer :: diag(:)
169 
170  ndof = hecmat%NDOF
171 
172  ! Extract diagonal components once (efficient, called only once)
173  diag => hecmw_mat_diag(hecmat)
174 
175  ! Loop over contact pairs
176  do i = 1, fstrsolid%n_contacts
177  grpid = fstrsolid%contacts(i)%group
178  if(.not. fstr_iscontactactive(fstrsolid, grpid, cstep)) cycle
179 
180  call calc_contact_pair_refstiff(fstrsolid%contacts(i), diag, ndof, hecmesh)
181  enddo
182 
183  ! Same for embeds
184  do i = 1, fstrsolid%n_embeds
185  grpid = fstrsolid%embeds(i)%group
186  if(.not. fstr_isembedactive(fstrsolid, grpid, cstep)) cycle
187 
188  call calc_contact_pair_refstiff(fstrsolid%embeds(i), diag, ndof, hecmesh)
189  enddo
190 
191  deallocate(diag)
192  end subroutine fstr_calc_contact_refstiff
193 
195  subroutine fstr_scan_contact_state( cstep, sub_step, cont_step, dt, ctAlgo, hecMESH, fstrSOLID, infoCTChange )
196  integer(kind=kint), intent(in) :: cstep
197  integer(kind=kint), intent(in) :: sub_step
198  integer(kind=kint), intent(in) :: cont_step
199  real(kind=kreal), intent(in) :: dt
200  integer(kind=kint), intent(in) :: ctAlgo
201  type( hecmwst_local_mesh ), intent(in) :: hecMESH
202  type(fstr_solid), intent(inout) :: fstrSOLID
203  type(fstr_info_contactchange), intent(inout):: infoCTChange
204  character(len=9) :: flag_ctAlgo
205  integer(kind=kint) :: i, grpid
206  logical :: iactive, is_init
207 
208  if( associated( fstrsolid%CONT_RELVEL ) ) fstrsolid%CONT_RELVEL(:) = 0.d0
209  if( associated( fstrsolid%CONT_STATE ) ) fstrsolid%CONT_STATE(:) = 0.d0
210 
211  if( ctalgo == kcaslagrange ) then
212  flag_ctalgo = 'SLagrange'
213  elseif( ctalgo == kcaalagrange ) then
214  flag_ctalgo = 'ALagrange'
215  endif
216 
217  ! P.A. We redefine fstrSOLID%ddunode as current coordinate of every nodes
218  ! fstrSOLID%ddunode(:) = fstrSOLID%unode(:) + fstrSOLID%dunode(:)
219  do i = 1, size(fstrsolid%unode)
220  fstrsolid%ddunode(i) = hecmesh%node(i) + fstrsolid%unode(i) + fstrsolid%dunode(i)
221  enddo
222  active = .false.
223 
224  infoctchange%contact2free = 0
225  infoctchange%contact2neighbor = 0
226  infoctchange%contact2diffLpos = 0
227  infoctchange%free2contact = 0
228  infoctchange%contactNode_current = 0
229 
230  is_init = ( cstep == 1 .and. sub_step == 1 .and. cont_step == 0 )
231 
232  do i=1,fstrsolid%n_contacts
233  grpid = fstrsolid%contacts(i)%group
234  if( .not. fstr_iscontactactive( fstrsolid, grpid, cstep ) ) then
235  call clear_contact_state(fstrsolid%contacts(i)); cycle
236  endif
237  call scan_contact_state( fstrsolid%contacts(i), &
238  & fstrsolid%ddunode(:), fstrsolid%dunode(:), infoctchange, &
239  & hecmesh%global_node_ID(:), hecmesh%global_elem_ID(:), &
240  & is_init, iactive, hecmesh=hecmesh, &
241  & flag_ctalgo=flag_ctalgo, ndforce=fstrsolid%QFORCE(:) )
242  if( .not. active ) active = iactive
243  enddo
244 
245  do i=1,fstrsolid%n_embeds
246  grpid = fstrsolid%embeds(i)%group
247  if( .not. fstr_isembedactive( fstrsolid, grpid, cstep ) ) then
248  call clear_contact_state(fstrsolid%embeds(i)); cycle
249  endif
250  call scan_embed_state( flag_ctalgo, fstrsolid%embeds(i), fstrsolid%ddunode(:), fstrsolid%dunode(:), &
251  & fstrsolid%QFORCE(:), infoctchange, hecmesh%global_node_ID(:), hecmesh%global_elem_ID(:), is_init, iactive )
252  if( .not. active ) active = iactive
253  enddo
254 
255  if( is_init .and. ctalgo == kcaslagrange .and. fstrsolid%n_contacts > 0 ) &
256  & call remove_duplication_tiedcontact( cstep, hecmesh, fstrsolid, infoctchange )
257 
258  infoctchange%contactNode_current = infoctchange%contactNode_previous+infoctchange%free2contact-infoctchange%contact2free
259  infoctchange%contactNode_previous = infoctchange%contactNode_current
260 
261  if( .not. active ) then
262  if( associated( fstrsolid%CONT_NFORCE ) ) fstrsolid%CONT_NFORCE(:) = 0.d0
263  if( associated( fstrsolid%CONT_FRIC ) ) fstrsolid%CONT_FRIC(:) = 0.d0
264  end if
265 
266  end subroutine
267 
269  subroutine fstr_scan_contact_state_exp( cstep, hecMESH, fstrSOLID, infoCTChange )
270  integer(kind=kint), intent(in) :: cstep
271  type( hecmwst_local_mesh ), intent(in) :: hecMESH
272  type(fstr_solid), intent(inout) :: fstrSOLID
273  type(fstr_info_contactchange), intent(inout) :: infoCTChange
274 
275  integer(kind=kint) :: i
276  logical :: iactive, is_init
277 
278 
279  ! P.A. We redefine fstrSOLID%ddunode as current coordinate of every nodes
280  ! fstrSOLID%ddunode(:) = fstrSOLID%unode(:) + fstrSOLID%dunode(:)
281  do i = 1, size(fstrsolid%unode)
282  fstrsolid%ddunode(i) = hecmesh%node(i) + fstrsolid%unode(i) + fstrsolid%dunode(i)
283  enddo
284  infoctchange%active = .false.
285 
286  infoctchange%contact2free = 0
287  infoctchange%contact2neighbor = 0
288  infoctchange%contact2diffLpos = 0
289  infoctchange%free2contact = 0
290  infoctchange%contactNode_current = 0
291 
292  is_init = ( cstep == 1 )
293 
294  do i=1,fstrsolid%n_contacts
295  ! grpid = fstrSOLID%contacts(i)%group
296  ! if( .not. fstr_isContactActive( fstrSOLID, grpid, cstep ) ) then
297  ! call clear_contact_state(fstrSOLID%contacts(i)); cycle
298  ! endif
299 
300  call scan_contact_state( fstrsolid%contacts(i), fstrsolid%ddunode(:), fstrsolid%dunode(:), &
301  & infoctchange, hecmesh%global_node_ID(:), hecmesh%global_elem_ID(:), is_init, iactive, hecmesh )
302 
303  infoctchange%active = infoctchange%active .or. iactive
304  enddo
305 
306  infoctchange%contactNode_current = infoctchange%contactNode_previous+infoctchange%free2contact-infoctchange%contact2free
307  infoctchange%contactNode_previous = infoctchange%contactNode_current
308  fstrsolid%ddunode = 0.d0
309  end subroutine
310 
311  logical function fstr_is_contact_active()
312  fstr_is_contact_active = active
313  end function
314 
315  subroutine fstr_set_contact_active( a )
316  logical, intent(in) :: a
317  active = a
318  end subroutine
319 
320  logical function fstr_is_contact_conv(ctAlgo,infoCTChange,hecMESH)
321  integer(kind=kint), intent(in) :: ctalgo
322  type (fstr_info_contactchange), intent(in) :: infoctchange
323  type (hecmwst_local_mesh), intent(in) :: hecmesh
324 
325  fstr_is_contact_conv = .false.
326  if( infoctchange%contact2free+infoctchange%contact2neighbor+ &
327  infoctchange%contact2difflpos+infoctchange%free2contact == 0 ) &
328  fstr_is_contact_conv = .true.
329 
330  call hecmw_allreduce_l1(hecmesh, fstr_is_contact_conv, hecmw_land)
331  end function
332 
333  logical function fstr_is_matrixstructure_changed(infoCTChange)
334  type (fstr_info_contactchange) :: infoctchange
336  if( infoctchange%contact2free+infoctchange%contact2neighbor+infoctchange%free2contact > 0 ) &
338  end function
339 
340  subroutine fstr_update_contact_multiplier( cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, ctchanged )
341  integer(kind=kint), intent(in) :: cstep
342  integer(kind=kint), intent(in) :: ctAlgo
343  type( hecmwst_local_mesh ), intent(in) :: hecMESH
344  type(hecmwst_matrix_lagrange), intent(in) :: hecLagMAT
345  type(fstr_solid), intent(inout) :: fstrSOLID
346  logical, intent(out) :: ctchanged
347 
348  integer(kind=kint) :: i, nc, algtype, grpid
349 
350  gnt = 0.d0; ctchanged = .false.
351  nc = fstrsolid%n_contacts+fstrsolid%n_embeds
352  do i=1, fstrsolid%n_contacts
353  grpid = fstrsolid%contacts(i)%group
354  if( .not. fstr_iscontactactive( fstrsolid, grpid, cstep ) ) cycle
355 
356  algtype = fstrsolid%contacts(i)%algtype
357  if( algtype == contactsslid .or. algtype == contactfslid ) then
358  call update_contact_multiplier( ctalgo, fstrsolid%contacts(i), hecmesh%node(:), fstrsolid%unode(:), &
359  fstrsolid%dunode(:), fstrsolid%contacts(i)%fcoeff, &
360  hecmesh, heclagmat, gnt, ctchanged )
361  else if( algtype == contacttied ) then
362  call update_tied_multiplier( fstrsolid%contacts(i), fstrsolid%unode(:), fstrsolid%dunode(:), &
363  & ctchanged )
364  endif
365  enddo
366 
367  do i=1, fstrsolid%n_embeds
368  call update_tied_multiplier( fstrsolid%embeds(i), fstrsolid%unode(:), fstrsolid%dunode(:), &
369  & ctchanged )
370  enddo
371 
372  call hecmw_allreduce_l1(hecmesh, ctchanged, hecmw_lor)
373 
374  if( nc>0 ) gnt = gnt/nc
375  end subroutine
376 
378  subroutine fstr_update_contact_tangentforce( cstep, fstrSOLID )
379  integer(kind=kint), intent(in) :: cstep
380  type(fstr_solid), intent(inout) :: fstrSOLID
381 
382  integer(kind=kint) :: i, grpid
383 
384  do i=1, fstrsolid%n_contacts
385 
386  grpid = fstrsolid%contacts(i)%group
387  if( .not. fstr_iscontactactive( fstrsolid, grpid, cstep ) ) cycle
388 
389  call update_contact_tangentforce( fstrsolid%contacts(i) )
390  enddo
391  end subroutine
392 
394  subroutine fstr_update_contact_state_vectors( fstrSOLID, dt )
395  type(fstr_solid), intent(inout) :: fstrsolid
396  real(kind=kreal), intent(in) :: dt
397 
398  integer(kind=kint) :: i
399 
400  if( associated( fstrsolid%CONT_RELVEL ) ) fstrsolid%CONT_RELVEL(:) = 0.d0
401  if( associated( fstrsolid%CONT_STATE ) ) fstrsolid%CONT_STATE(:) = 0.d0
402 
403  do i=1, fstrsolid%n_contacts
404  call update_contact_state_vectors( fstrsolid%contacts(i), dt, fstrsolid%CONT_RELVEL, fstrsolid%CONT_STATE )
405  enddo
406  end subroutine
407 
408 end module mcontact
real(kind=kreal) function, dimension(:), pointer, public hecmw_mat_diag(hecMAT)
Extract diagonal components from matrix D into a 1D vector Returns: diag(i) = D(ndof*ndof*(node-1) + ...
Definition: hecmw.f90:6
Contact processing at assembly level (all pairs in one tContact object)
subroutine calcu_contact_ndforce_nodesurf(purpose, ctAlgo, contact, coord, disp, ddisp, lagrange_array, conMAT, CONT_NFORCE, CONT_FRIC, hecLagMAT)
This subroutine calculates contact nodal force for each contact pair and assembles it into contact ma...
subroutine, public calc_contact_pair_refstiff(contact, diag, ndof, hecMESH)
Calculate reference stiffness for one contact pair.
subroutine update_contact_multiplier(ctAlgo, contact, coord, disp, ddisp, fcoeff, hecMESH, hecLagMAT, gnt, ctchanged)
This subroutine update lagrangian multiplier and the distance between contacting nodes.
subroutine calcu_contact_stiffness_nodesurf(ctAlgo, contact, coord, disp, ddisp, iter, lagrange_array, conMAT, hecLagMAT)
This subroutine calculates contact stiffness for each contact pair and assembles it into global stiff...
subroutine update_contact_tangentforce(contact)
subroutine update_tied_multiplier(contact, disp, ddisp, ctchanged)
This subroutine update lagrangian multiplier and the distance between contacting nodes.
Contact mechanics calculations at element level (single contact pair)
This module provides geometric calculations for contact.
MPC (Multi-Point Constraint) processing for contact analysis.
Contact output vector processing (initialization, area calculation, parallel)
subroutine, public fstr_setup_parancon_contactvalue(hecMESH, ndof, vec, vtype)
subroutine, public calc_contact_area(hecMESH, fstrSOLID, flag)
subroutine, public initialize_contact_output_vectors(fstrSOLID, hecMAT)
subroutine, public setup_contact_elesurf_for_area(cstep, hecMESH, fstrSOLID)
subroutine, public initialize_embed_vectors(fstrSOLID, hecMAT)
subroutine, public update_contact_state_vectors(contact, dt, relvel_vec, state_vec)
Update contact state output vectors (CONT_RELVEL, CONT_STATE)
Contact search and state scanning at search level (all pairs in one tContact)
subroutine scan_contact_state(contact, currpos, currdisp, infoCTChange, nodeID, elemID, is_init, active, hecMESH, flag_ctAlgo, ndforce)
This subroutine update contact states, which include.
subroutine remove_duplication_tiedcontact(cstep, hecMESH, fstrSOLID, infoCTChange)
Scanning contact state.
subroutine scan_embed_state(flag_ctAlgo, embed, currpos, currdisp, ndforce, infoCTChange, nodeID, elemID, is_init, active, B)
This subroutine update contact states, which include.
This module defines common data and basic structures for analysis.
Definition: m_fstr.F90:15
logical function fstr_isembedactive(fstrSOLID, nbc, cstep)
Definition: m_fstr.F90:1082
logical function fstr_iscontactactive(fstrSOLID, nbc, cstep)
Definition: m_fstr.F90:1072
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
Definition: m_fstr.F90:61
real(kind=kreal) dt
ANALYSIS CONTROL for NLGEOM and HEAT.
Definition: m_fstr.F90:141
integer(kind=kint), parameter kcaalagrange
Definition: m_fstr.F90:62
Top-level contact analysis module (System level)
subroutine fstr_set_contact_active(a)
subroutine fstr_update_ndforce_contact(cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT)
Compute contact forces for residual vector (conMATB).
logical function fstr_is_matrixstructure_changed(infoCTChange)
logical function fstr_is_contact_conv(ctAlgo, infoCTChange, hecMESH)
subroutine fstr_contact_ndforce_core(purpose, cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT)
Core routine: compute contact nodal forces for all contact/embed pairs. purpose == kctForResidual: as...
subroutine fstr_addcontactstiffness(cstep, ctAlgo, iter, hecMESH, conMAT, hecLagMAT, fstrSOLID)
subroutine fstr_calc_contact_output_force(cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT)
Compute contact forces for output (CONT_NFORCE/CONT_FRIC).
subroutine fstr_update_contact_multiplier(cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, ctchanged)
subroutine, public fstr_update_contact_state_vectors(fstrSOLID, dt)
Update contact state output vectors for all contacts.
subroutine fstr_calc_contact_refstiff(cstep, hecMESH, hecMAT, fstrSOLID)
Calculate reference stiffness for all contact pairs (System Level)
subroutine fstr_scan_contact_state_exp(cstep, hecMESH, fstrSOLID, infoCTChange)
Scanning contact state.
subroutine fstr_scan_contact_state(cstep, sub_step, cont_step, dt, ctAlgo, hecMESH, fstrSOLID, infoCTChange)
Scanning contact state.
logical function fstr_is_contact_active()
subroutine fstr_update_contact_tangentforce(cstep, fstrSOLID)
Update tangent force.
This module manages the data structure for contact calculation.
integer, parameter contactsslid
integer, parameter kctforresidual
purpose flag for contact force calculation
real(kind=kreal), dimension(2), save gnt
1:current average penetration; 2:current relative tangent displacement
subroutine clear_contact_state(contact)
Reset contact state all to free.
integer, parameter contacttied
contact type or algorithm definition
integer, parameter kctforoutput
compute contact force for output (CONT_NFORCE/CONT_FRIC)
integer, parameter contactfslid