FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_contact_assembly.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 !-------------------------------------------------------------------------------
8  use hecmw
9  use m_fstr
10  use mcontactdef
14  implicit none
15 
17 
18 contains
19 
21  subroutine calc_contact_pair_refstiff(contact, diag, ndof, hecMESH)
22  type(tcontact), intent(inout) :: contact
23  real(kind=kreal), intent(in) :: diag(:)
24  integer(kind=kint), intent(in) :: ndof
25  type(hecmwst_local_mesh), intent(in) :: hecmesh
26 
27  integer(kind=kint) :: j, slave_node, master_node, nnode, ctsurf
28  integer(kind=kint) :: idx_start, idx_end
29  real(kind=kreal) :: maxv
30 
31  maxv = 0.0d0
32 
33  ! Loop over slave nodes
34  do j = 1, size(contact%slave)
35  slave_node = contact%slave(j)
36  idx_start = ndof * (slave_node - 1) + 1
37  idx_end = ndof * slave_node
38  maxv = max(maxv, maxval(diag(idx_start:idx_end)))
39  enddo
40 
41  ! Loop over master surfaces and nodes
42  do ctsurf = 1, size(contact%master)
43  nnode = size(contact%master(ctsurf)%nodes)
44  do j = 1, nnode
45  master_node = contact%master(ctsurf)%nodes(j)
46  idx_start = ndof * (master_node - 1) + 1
47  idx_end = ndof * master_node
48  maxv = max(maxv, maxval(diag(idx_start:idx_end)))
49  enddo
50  enddo
51 
52  ! Parallel reduction
53  call hecmw_allreduce_r1(hecmesh, maxv, hecmw_max)
54 
55  ! Set reference stiffness for this contact pair
56  contact%refStiff = maxv
57 
58  ! Report penalty settings
59  if (hecmw_comm_get_rank() == 0) then
60  write(*,'(A,A,A,1pE12.3,A,1pE12.3,A,1pE12.3)') " Contact [", &
61  trim(contact%pair_name), "] set penalty: normal & tied ", &
62  contact%nPenalty * contact%refStiff, ", tangential ", &
63  contact%tPenalty * contact%refStiff, ", refStiff ", contact%refStiff
64  endif
65 
66  end subroutine calc_contact_pair_refstiff
67 
69  subroutine assemble_contact_force_residual(nnode,ndLocal,id_lagrange,ctNForce,ctTForce,conMAT)
70  integer(kind=kint), intent(in) :: nnode
71  integer(kind=kint), intent(in) :: ndLocal(nnode + 1)
72  integer(kind=kint), intent(in) :: id_lagrange
73  real(kind=kreal), intent(in) :: ctnforce((nnode+1)*3+1)
74  real(kind=kreal), intent(in) :: cttforce((nnode+1)*3+1)
75  type(hecmwst_matrix), intent(inout) :: conmat
76 
77  integer(kind=kint) :: i, inod, idx
78 
79  do i = 1, nnode + 1
80  inod = ndlocal(i)
81  idx = (inod-1)*3+1
82  conmat%B(idx:idx+2) = conmat%B(idx:idx+2) + ctnforce((i-1)*3+1:(i-1)*3+3) + cttforce((i-1)*3+1:(i-1)*3+3)
83  enddo
84 
85  if( id_lagrange > 0 ) then
86  conmat%B(conmat%NP*conmat%NDOF+id_lagrange) = ctnforce((nnode+1)*3+1) + cttforce((nnode+1)*3+1)
87  endif
88 
89  end subroutine assemble_contact_force_residual
90 
92  subroutine assemble_contact_force_output(nnode,ndLocal,ctNForce,ctTForce,cont_nforce,cont_fric)
93  integer(kind=kint), intent(in) :: nnode
94  integer(kind=kint), intent(in) :: ndlocal(nnode + 1)
95  real(kind=kreal), intent(in) :: ctnforce((nnode+1)*3+1)
96  real(kind=kreal), intent(in) :: cttforce((nnode+1)*3+1)
97  real(kind=kreal), pointer, intent(inout) :: cont_nforce(:)
98  real(kind=kreal), pointer, optional, intent(inout) :: cont_fric(:)
99 
100  integer(kind=kint) :: i, inod, idx
101 
102  do i = 1, nnode + 1
103  inod = ndlocal(i)
104  idx = (inod-1)*3+1
105  cont_nforce(idx:idx+2) = cont_nforce(idx:idx+2) + ctnforce((i-1)*3+1:(i-1)*3+3)
106  if( present(cont_fric) ) cont_fric(idx:idx+2) = cont_fric(idx:idx+2) + cttforce((i-1)*3+1:(i-1)*3+3)
107  enddo
108 
109  end subroutine assemble_contact_force_output
110 
113  subroutine update_contact_multiplier( ctAlgo, contact, coord, disp, ddisp, fcoeff, &
114  hecMESH, hecLagMAT, gnt, ctchanged )
115  integer(kind=kint), intent(in) :: ctAlgo
116  type( tcontact ), intent(inout) :: contact
117  real(kind=kreal), intent(in) :: coord(:)
118  real(kind=kreal), intent(in) :: disp(:)
119  real(kind=kreal), intent(in) :: ddisp(:)
120  real(kind=kreal), intent(in) :: fcoeff
121  type(hecmwst_local_mesh), intent(in) :: hecmesh
122  type(hecmwst_matrix_lagrange), intent(in) :: hecLagMAT
123  real(kind=kreal), intent(out) :: gnt(2)
124  logical, intent(inout) :: ctchanged
125 
126  integer(kind=kint) :: slave, etype, master
127  integer(kind=kint) :: nn, i, cnt
128  real(kind=kreal) :: lgnt(2)
129  integer(kind=kint) :: ndLocal(l_max_elem_node+1)
130  real(kind=kreal) :: ctnforce(l_max_elem_node*3+3)
131  real(kind=kreal) :: cttforce(l_max_elem_node*3+3)
132  real(kind=kreal) :: max_jump_ratio, jump_ratio_local
133  real(kind=kreal) :: mut_old, mut_new, threthold
134 
135  cnt = 0
136  lgnt(:) = 0.d0
137  max_jump_ratio = 0.0d0
138 
139  do i = 1, size(contact%slave)
140  if(contact%states(i)%state == contactfree) cycle ! not in contact
141 
142  slave = contact%slave(i)
143  master = contact%states(i)%surface
144  nn = size(contact%master(master)%nodes)
145  etype = contact%master(master)%etype
146 
147  ndlocal(1) = slave
148  ndlocal(2:nn+1) = contact%master(master)%nodes(1:nn)
149 
150  ! Update multiplier and calculate forces
151  call updatecontactmultiplier_alag(contact%states(i), ndlocal(1:nn+1), coord, disp, ddisp, &
152  contact%nPenalty * contact%refStiff, contact%tPenalty * contact%refStiff, &
153  fcoeff, etype, lgnt, ctchanged, ctnforce, cttforce, jump_ratio_local)
154 
155  ! Track maximum jump ratio
156  max_jump_ratio = max(max_jump_ratio, jump_ratio_local)
157 
158  cnt = cnt + 1
159  enddo
160 
161  if(cnt > 0) lgnt(:) = lgnt(:) / cnt
162  gnt = gnt + lgnt
163 
164  call hecmw_allreduce_r1(hecmesh, max_jump_ratio, hecmw_max)
165 
166  ! Adjust tPenalty
167  threthold = 100.d0
168  if (max_jump_ratio > threthold) then
169  mut_old = contact%tPenalty * contact%refStiff
170  contact%tPenalty = contact%tPenalty * max(1.d0/dsqrt(threthold), 1.0d0/dsqrt(max_jump_ratio))
171  mut_new = contact%tPenalty * contact%refStiff
172  if (hecmw_comm_get_rank() == 0) then
173  write(*,'(A,A,A,1pE12.3,A,1pE12.3,A)') " Contact [", trim(contact%pair_name), &
174  "] tangential penalty adjusted: ", mut_old, " -> ", mut_new, " (friction jump)"
175  endif
176  endif
177 
178  end subroutine update_contact_multiplier
179 
182  subroutine update_tied_multiplier( contact, disp, ddisp, ctchanged )
183  type( tcontact ), intent(inout) :: contact
184  real(kind=kreal), intent(in) :: disp(:)
185  real(kind=kreal), intent(in) :: ddisp(:)
186  logical, intent(inout) :: ctchanged
187 
188  integer(kind=kint) :: slave, etype, master
189  integer(kind=kint) :: nn, i, j, iss
190  real(kind=kreal) :: dg(3), dgmax
191  real(kind=kreal) :: shapefunc(l_max_surface_node)
192  real(kind=kreal) :: edisp(3*l_max_elem_node+3)
193  real(kind=kreal) :: mu
194 
195  ! Calculate penalty from contact structure
196  mu = contact%nPenalty * contact%refStiff
197 
198  do i= 1, size(contact%slave)
199  if( contact%states(i)%state==contactfree ) cycle ! not in contact
200  slave = contact%slave(i)
201  edisp(1:3) = disp(3*slave-2:3*slave)+ddisp(3*slave-2:3*slave)
202  master = contact%states(i)%surface
203 
204  nn = size( contact%master(master)%nodes )
205  etype = contact%master(master)%etype
206  do j=1,nn
207  iss = contact%master(master)%nodes(j)
208  edisp(3*j+1:3*j+3) = disp(3*iss-2:3*iss)+ddisp(3*iss-2:3*iss)
209  enddo
210  call getshapefunc( etype, contact%states(i)%lpos(1:2), shapefunc )
211 
212  ! normal component
213  dg(1:3) = edisp(1:3)
214  do j=1,nn
215  dg(1:3) = dg(1:3)-shapefunc(j)*edisp(3*j+1:3*j+3)
216  enddo
217 
218  contact%states(i)%multiplier(1:3) = contact%states(i)%multiplier(1:3) + mu*dg(1:3)
219 
220  ! check if tied constraint converged
221  dgmax = 0.d0
222  do j=1,(nn+1)*3
223  dgmax = dgmax + dabs(edisp(j))
224  enddo
225  dgmax = dgmax/dble((nn+1)*3)
226  do j=1,3
227  if( dabs(dg(j))/dmax1(1.d0,dgmax) > 1.d-3 ) ctchanged = .true.
228  enddo
229 
230  enddo
231  end subroutine
232 
233  subroutine update_contact_tangentforce( contact )
234  type( tcontact ), intent(inout) :: contact
235 
236  integer(kind=kint) :: i
237 
238  do i= 1, size(contact%slave)
239  if( contact%states(i)%state==contactfree ) then
240  contact%states(i)%tangentForce(1:3) = 0.d0
241  contact%states(i)%tangentForce_trial(1:3) = 0.d0
242  contact%states(i)%tangentForce_final(1:3) = 0.d0
243  else
244  contact%states(i)%tangentForce(1:3) = contact%states(i)%tangentForce_final(1:3)
245  end if
246  contact%states(i)%tangentForce1(1:3) = contact%states(i)%tangentForce(1:3)
247  enddo
248  end subroutine update_contact_tangentforce
249 
252  subroutine calcu_contact_stiffness_nodesurf( ctAlgo, contact, coord, disp, iter, lagrange_array, &
253  conMAT, hecLagMAT)
254  integer(kind=kint), intent(in) :: ctalgo
255  type(tcontact), intent(inout) :: contact
256  real(kind=kreal), intent(in) :: coord(:)
257  real(kind=kreal), intent(in) :: disp(:)
258  integer(kind=kint), intent(in) :: iter
259  real(kind=kreal), intent(in) :: lagrange_array(:)
260  type(hecmwst_matrix), intent(inout) :: conmat
261  type(hecmwst_matrix_lagrange), intent(inout) :: hecLagMAT
262 
263  integer(kind=kint) :: ctsurf, nnode, ndLocal(21), etype
264  integer(kind=kint) :: j, k, algtype, id_lagrange
265  real(kind=kreal) :: lagrange
266  real(kind=kreal) :: stiffness((l_max_surface_node+1)*3+1, (l_max_surface_node+1)*3+1)
267  real(kind=kreal) :: elecoord(3, l_max_surface_node)
268  real(kind=kreal) :: force(l_max_surface_node*3+3)
269 
270  algtype = contact%algtype
271 
272  do j = 1, size(contact%slave)
273 
274  if( contact%states(j)%state == contactfree ) cycle
275 
276  ctsurf = contact%states(j)%surface
277  etype = contact%master(ctsurf)%etype
278  nnode = size(contact%master(ctsurf)%nodes)
279  ndlocal(1) = contact%slave(j)
280  ndlocal(2:nnode+1) = contact%master(ctsurf)%nodes(1:nnode)
281 
282  ! Prepare master node coordinates for ALagrange (deformed configuration)
283  do k = 1, nnode
284  elecoord(1:3, k) = coord(3*ndlocal(k+1)-2:3*ndlocal(k+1)) + disp(3*ndlocal(k+1)-2:3*ndlocal(k+1))
285  enddo
286 
287  if( algtype == contactsslid .or. algtype == contactfslid ) then
288 
289  if( ctalgo == kcaslagrange ) then
290  id_lagrange = heclagmat%lag_node_table(ndlocal(1)) - 1
291  id_lagrange = id_lagrange + 1
292  lagrange = lagrange_array(id_lagrange)
293  call getcontactstiffness_slag(contact%states(j), contact%master(ctsurf), iter, &
294  contact%tPenalty, contact%fcoeff, lagrange, stiffness)
295 
296  ! Assemble contact stiffness matrix of contact pair into global stiffness matrix
297  call hecmw_mat_ass_contactlag(nnode, ndlocal, id_lagrange, contact%fcoeff, stiffness, conmat, heclagmat)
298 
299  else if( ctalgo == kcaalagrange ) then
300  call getcontactstiffness_alag(contact%states(j), contact%master(ctsurf), elecoord(:,1:nnode), &
301  contact%nPenalty * contact%refStiff, contact%tPenalty * contact%refStiff, &
302  contact%fcoeff, contact%symmetric, stiffness, force)
303 
304  ! Assemble contact stiffness matrix into global stiffness matrix
305  call hecmw_mat_ass_elem(conmat, nnode+1, ndlocal, stiffness)
306 
307  end if
308 
309 
310  else if( algtype == contacttied ) then
311 
312  if( ctalgo == kcaslagrange ) then
313  id_lagrange = heclagmat%lag_node_table(ndlocal(1)) - 1
314  do k = 1, 3
315  id_lagrange = id_lagrange + 1
316  lagrange = lagrange_array(id_lagrange)
317 
318  call gettiedstiffness_slag(contact%states(j), contact%master(ctsurf), k, stiffness)
319  ! Assemble contact stiffness matrix of contact pair into global stiffness matrix
320  call hecmw_mat_ass_contactlag(nnode, ndlocal, id_lagrange, 0.d0, stiffness, conmat, heclagmat)
321  enddo
322 
323  else if( ctalgo == kcaalagrange ) then
324  call gettiedstiffness_alag(contact%states(j), contact%master(ctsurf), &
325  contact%nPenalty * contact%refStiff, stiffness, force)
326 
327  ! Assemble contact stiffness matrix into global stiffness matrix
328  call hecmw_mat_ass_elem(conmat, nnode+1, ndlocal, stiffness)
329 
330  end if
331 
332  endif
333 
334  enddo
335 
336  end subroutine calcu_contact_stiffness_nodesurf
337 
342  subroutine calcu_contact_ndforce_nodesurf( purpose, ctAlgo, contact, coord, disp, ddisp, lagrange_array, &
343  conMAT, CONT_NFORCE, CONT_FRIC, hecLagMAT )
344  integer(kind=kint), intent(in) :: purpose
345  integer(kind=kint), intent(in) :: ctAlgo
346  type( tcontact ), intent(inout) :: contact
347  real(kind=kreal), intent(in) :: coord(:)
348  real(kind=kreal), intent(in) :: disp(:)
349  real(kind=kreal), intent(in) :: ddisp(:)
350  real(kind=kreal), intent(in) :: lagrange_array(:)
351  type(hecmwst_matrix), intent(inout) :: conmat
352  real(kind=kreal), pointer :: cont_nforce(:)
353  real(kind=kreal), pointer :: cont_fric(:)
354  type(hecmwst_matrix_lagrange), intent(in) :: heclagmat
355 
356  integer(kind=kint) :: ctsurf, nnode, ndlocal(21)
357  integer(kind=kint) :: j, k, algtype, id_lagrange
358  real(kind=kreal) :: ndcoord(21*3)
359  real(kind=kreal) :: ndu(21*3), nddu(21*3)
360  real(kind=kreal) :: lagrange
361  real(kind=kreal) :: ctnforce(21*3+1)
362  real(kind=kreal) :: cttforce(21*3+1)
363  real(kind=kreal) :: mu_n, mu_t
364  logical :: if_flag
365  real(kind=kreal) :: ctime, etime
366  integer(kind=kint) :: if_type
367 
368  algtype = contact%algtype
369  if_flag = (contact%if_type /= 0)
370  if(if_flag)then
371  ctime = contact%ctime
372  etime = contact%if_etime
373  if_type = contact%if_type
374  end if
375 
376  do j = 1, size(contact%slave)
377 
378  if( contact%states(j)%state == contactfree ) cycle
379  if(if_flag) call set_shrink_factor(ctime, contact%states(j), etime, if_type)
380 
381  ctsurf = contact%states(j)%surface
382  nnode = size(contact%master(ctsurf)%nodes)
383  ndlocal(1) = contact%slave(j)
384  ndlocal(2:nnode+1) = contact%master(ctsurf)%nodes(1:nnode)
385  do k = 1, nnode+1
386  nddu((k-1)*3+1:(k-1)*3+3) = ddisp((ndlocal(k)-1)*3+1:(ndlocal(k)-1)*3+3)
387  ndu((k-1)*3+1:(k-1)*3+3) = disp((ndlocal(k)-1)*3+1:(ndlocal(k)-1)*3+3) + nddu((k-1)*3+1:(k-1)*3+3)
388  ndcoord((k-1)*3+1:(k-1)*3+3) = coord((ndlocal(k)-1)*3+1:(ndlocal(k)-1)*3+3) + ndu((k-1)*3+1:(k-1)*3+3)
389  enddo
390 
391  ! --- Determine penalty parameters: zero for output (multiplier-only)
392  if( ctalgo == kcaalagrange .and. purpose == kctforoutput ) then
393  mu_n = 0.0d0
394  mu_t = 0.0d0
395  else
396  mu_n = contact%nPenalty * contact%refStiff
397  mu_t = contact%tPenalty * contact%refStiff
398  endif
399 
400  if( algtype == contactsslid .or. algtype == contactfslid ) then
401  ! Obtain contact nodal force vector of contact pair
402  if(if_flag) call get_shrink_elemact_surf(contact%states(j),ndcoord, nnode)
403 
404  if( ctalgo == kcaslagrange ) then
405  id_lagrange = heclagmat%lag_node_table(ndlocal(1)) - 1
406  id_lagrange = id_lagrange + 1
407  lagrange = lagrange_array(id_lagrange)
408  call getcontactnodalforce_slag(contact%states(j),contact%master(ctsurf),ndcoord,nddu, &
409  contact%tPenalty,contact%fcoeff,lagrange,ctnforce,cttforce,.true.)
410 
411  else if( ctalgo == kcaalagrange ) then
412  id_lagrange = 0
413  lagrange = 0.d0
414  call getcontactnodalforce_alag(contact%states(j),contact%master(ctsurf),ndcoord,nddu, &
415  mu_n, mu_t, contact%fcoeff,lagrange,ctnforce,cttforce,.true.)
416 
417  end if
418 
419  ! Assemble contact force
420  if( purpose == kctforresidual ) then
421  call assemble_contact_force_residual(nnode,ndlocal,id_lagrange,ctnforce,cttforce,conmat)
422  else
423  call assemble_contact_force_output(nnode,ndlocal,ctnforce,cttforce,cont_nforce,cont_fric)
424  endif
425 
426  else if( algtype == contacttied ) then
427 
428  if( ctalgo == kcaslagrange ) then
429  id_lagrange = heclagmat%lag_node_table(ndlocal(1)) - 1
430  do k=1,3
431  id_lagrange = id_lagrange + 1
432  lagrange = lagrange_array(id_lagrange)
433  contact%states(j)%multiplier(k) = lagrange
434 
435  call gettiednodalforce_slag(contact%states(j),contact%master(ctsurf),k,ndu, &
436  & lagrange,ctnforce,cttforce)
437  if( purpose == kctforresidual ) then
438  call assemble_contact_force_residual(nnode,ndlocal,id_lagrange,ctnforce,cttforce,conmat)
439  else
440  call assemble_contact_force_output(nnode,ndlocal,ctnforce,cttforce,cont_nforce)
441  endif
442  end do
443 
444  else if( ctalgo == kcaalagrange ) then
445  id_lagrange = 0
446  call gettiednodalforce_alag(contact%states(j),contact%master(ctsurf),ndu, &
447  mu_n, ctnforce,cttforce)
448  if( purpose == kctforresidual ) then
449  call assemble_contact_force_residual(nnode,ndlocal,id_lagrange,ctnforce,cttforce,conmat)
450  else
451  call assemble_contact_force_output(nnode,ndlocal,ctnforce,cttforce,cont_nforce)
452  endif
453 
454  end if
455 
456  endif
457 
458  enddo
459 
460  end subroutine calcu_contact_ndforce_nodesurf
461 
462 end module m_fstr_contact_assembly
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 calcu_contact_stiffness_nodesurf(ctAlgo, contact, coord, disp, iter, lagrange_array, conMAT, hecLagMAT)
This subroutine calculates contact stiffness for each contact pair and assembles it into global stiff...
subroutine, public calc_contact_pair_refstiff(contact, diag, ndof, hecMESH)
Calculate reference stiffness for one contact pair.
subroutine assemble_contact_force_residual(nnode, ndLocal, id_lagrange, ctNForce, ctTForce, conMAT)
Assemble contact nodal force into residual vector (conMATB).
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 assemble_contact_force_output(nnode, ndLocal, ctNForce, ctTForce, cont_nforce, cont_fric)
Accumulate contact nodal force into output arrays (CONT_NFORCE/CONT_FRIC).
subroutine update_contact_tangentforce(contact)
subroutine update_tied_multiplier(contact, disp, ddisp, ctchanged)
This subroutine update lagrangian multiplier and the distance between contacting nodes.
Alag method implementations for contact element calculations.
subroutine, public updatecontactmultiplier_alag(ctState, ndLocal, coord, disp, ddisp, mu, mut, fcoeff, etype, lgnt, ctchanged, ctNForce, ctTForce, jump_ratio)
subroutine, public gettiedstiffness_alag(cstate, tSurf, mu, stiff, force)
subroutine, public getcontactnodalforce_alag(ctState, tSurf, ndCoord, ndDu, mu, mut, fcoeff, lagrange, ctNForce, ctTForce, cflag)
subroutine, public gettiednodalforce_alag(ctState, tSurf, ndu, mu, ctNForce, ctTForce)
subroutine, public getcontactstiffness_alag(cstate, tSurf, ele, mu, mut, fcoeff, symm, stiff, force)
Contact mechanics calculations at element level (single contact pair)
This module provides interference fit (shrink) functions for contact.
subroutine get_shrink_elemact_surf(cstate, coords, nnode)
subroutine, public set_shrink_factor(ctime, cstate, etime, if_type)
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
real(kind=kreal) etime
Definition: m_fstr.f90:141
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
Definition: m_fstr.f90:60
integer(kind=kint), parameter kcaalagrange
Definition: m_fstr.f90:61
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
integer, parameter contacttied
contact type or algorithm definition
integer, parameter kctforoutput
compute contact force for output (CONT_NFORCE/CONT_FRIC)
integer, parameter contactfree
contact state definition
integer, parameter contactfslid
Structure to includes all info needed by contact calculation.