FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_mat_con_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 !-------------------------------------------------------------------------------
6 ! stiffness matrix structure for the contact analysis
7 ! employing standard Lagrange multiplier algorithm
8 
10 
11  use m_fstr
12  use elementinfo
14 
15  implicit none
16  private
18  public :: hecmwst_matrix_lagrange
20  public :: fstr_mat_con_contact
22 
23  integer(kind=kint), save :: NPL_org, NPU_org
24  type(nodeRelated), pointer, save :: list_nodeRelated_org(:) => null()
25 
26  type(nodeRelated), pointer :: list_nodeRelated(:) => null()
27 
28  logical :: permission = .false.
29 
30 contains
31 
32  integer(kind=kint) function fstr_get_num_lagrange_pernode(algtype)
33  integer(kind=kint) :: algtype
34  if( algtype == contactsslid .or. algtype == contactfslid ) then
36  else if( algtype == contacttied ) then
38  endif
39  end function
40 
43 
44  type(hecmwst_matrix) :: hecmat
45 
46  if( associated(list_noderelated_org) ) return
47  call hecmw_construct_noderelated_from_hecmat(hecmat, npl_org, npu_org, list_noderelated_org)
48 
50 
53  subroutine fstr_mat_con_contact(cstep,contact_algo,hecMAT,fstrSOLID,hecLagMAT,infoCTChange,conMAT,is_contact_active_flag)
54 
55  integer(kind=kint) :: cstep
56  integer(kind=kint) :: contact_algo
57  type(hecmwst_matrix) :: hecmat
58  type(fstr_solid) :: fstrsolid
59  type(hecmwst_matrix_lagrange) :: heclagmat
60  type(fstr_info_contactchange) :: infoctchange
61 
62  integer(kind=kint) :: num_lagrange
63  integer(kind=kint) :: countnon0lu_node, countnon0lu_lagrange
64  integer(kind=kint) :: numnon0_node, numnon0_lagrange
66  type (hecmwst_matrix) :: conmat
67  logical, intent(in) :: is_contact_active_flag
68 
69  integer(kind=kint) :: i, j, grpid
70  integer(kind=kint) :: nlag
71 
72  num_lagrange = 0
73  if( contact_algo == kcaslagrange ) then
74  do i = 1, fstrsolid%n_contacts
75  grpid = fstrsolid%contacts(i)%group
76  if( .not. fstr_iscontactactive( fstrsolid, grpid, cstep ) ) cycle
77  nlag = fstr_get_num_lagrange_pernode(fstrsolid%contacts(i)%algtype)
78  do j = 1, size(fstrsolid%contacts(i)%slave)
79  if( .not. is_contact_active(fstrsolid%contacts(i)%states(j)%state) ) cycle
80  num_lagrange = num_lagrange + nlag
81  enddo
82  enddo
83 
84  do i = 1, fstrsolid%n_embeds
85  grpid = fstrsolid%embeds(i)%group
86  if( .not. fstr_iscontactactive( fstrsolid, grpid, cstep ) ) cycle
87  nlag = 3
88  do j = 1, size(fstrsolid%embeds(i)%slave)
89  if( .not. is_contact_active(fstrsolid%embeds(i)%states(j)%state) ) cycle
90  num_lagrange = num_lagrange + nlag
91  enddo
92  enddo
93  endif
94 
95  ! Get original list of related nodes
96  call hecmw_init_noderelated_from_org(hecmat%NP,num_lagrange,is_contact_active_flag,list_noderelated_org,list_noderelated)
97 
98  ! Construct new list of related nodes and Lagrange multipliers
99  countnon0lu_node = npl_org + npu_org
100  countnon0lu_lagrange = 0
101  if( is_contact_active_flag ) call getnewlistofrelatednodesandlagrangemultipliers(cstep,contact_algo, &
102  & hecmat%NP,fstrsolid,countnon0lu_node,countnon0lu_lagrange,list_noderelated)
103 
104  ! Construct new matrix structure(hecMAT&hecLagMAT)
105  numnon0_node = countnon0lu_node/2
106  numnon0_lagrange = countnon0lu_lagrange/2
107  call hecmw_construct_hecmat_from_noderelated(hecmat%N, hecmat%NP, hecmat%NDOF, &
108  & numnon0_node, num_lagrange, list_noderelated, hecmat)
109  call hecmw_construct_hecmat_from_noderelated(hecmat%N, hecmat%NP, hecmat%NDOF, &
110  & numnon0_node, num_lagrange, list_noderelated, conmat)
111  if( contact_algo == kcaslagrange ) call hecmw_construct_heclagmat_from_noderelated(hecmat%NP, &
112  & hecmat%NDOF, num_lagrange, numnon0_lagrange, is_contact_active_flag, list_noderelated, heclagmat)
113  call hecmw_finalize_noderelated(list_noderelated)
114 
115  ! Copy Lagrange multipliers
116  if( is_contact_active_flag .and. contact_algo == kcaslagrange ) &
117  call fstr_copy_lagrange_contact(fstrsolid,heclagmat)
118 
119  end subroutine fstr_mat_con_contact
120 
122  subroutine getnewlistofrelatednodesandlagrangemultipliers( &
123  & cstep, contact_algo, np, fstrSOLID, countNon0LU_node, countNon0LU_lagrange, list_nodeRelated )
124  integer(kind=kint),intent(in) :: cstep
125  integer(kind=kint),intent(in) :: contact_algo
126  integer(kind=kint),intent(in) :: np
127  type(fstr_solid),intent(in) :: fstrsolid
128  integer(kind=kint), intent(inout) :: countnon0lu_node, countnon0lu_lagrange
129  type(noderelated), pointer, intent(inout) :: list_noderelated(:)
130 
131  integer(kind=kint) :: grpid
132  integer(kind=kint) :: count_lagrange
133  integer(kind=kint) :: ctsurf, etype, nnode, ndlocal(l_max_surface_node + 1)
134  integer(kind=kint) :: i, j, k, nlag, algtype
135  real(kind=kreal) :: fcoeff
136  logical :: necessary_to_insert_node, necessary_to_insert_node_pair
137  logical :: is_contact_active_flag, is_damping_active_flag
138 
139  count_lagrange = 0
140  do i = 1, fstrsolid%n_contacts
141 
142  grpid = fstrsolid%contacts(i)%group
143  if( .not. fstr_iscontactactive( fstrsolid, grpid, cstep ) ) cycle
144 
145  fcoeff = fstrsolid%contacts(i)%fcoeff
146  necessary_to_insert_node = ( fcoeff /= 0.0d0 .or. contact_algo == kcaalagrange )
147 
148  algtype = fstrsolid%contacts(i)%algtype
149  nlag = fstr_get_num_lagrange_pernode(algtype)
150  if( contact_algo == kcaalagrange ) nlag = 1
151  if( algtype == contacttied ) permission = .true.
152 
153  do j = 1, size(fstrsolid%contacts(i)%slave)
154  ! stick or sliding contact is active
155  is_contact_active_flag = is_contact_active(fstrsolid%contacts(i)%states(j)%state)
156  ! damping is active
157  is_damping_active_flag = fstrsolid%contacts(i)%states(j)%state == contactnear .and. &
158  & is_damping_enabled(fstrsolid%contacts(i))
159 
160  if( is_contact_active_flag .or. is_damping_active_flag ) then
161 
162  ctsurf = fstrsolid%contacts(i)%states(j)%surface
163  etype = fstrsolid%contacts(i)%master(ctsurf)%etype
164  if( etype/=fe_tri3n .and. etype/=fe_quad4n ) &
165  stop " ##Error: This element type is not supported in contact analysis !!! "
166  nnode = size(fstrsolid%contacts(i)%master(ctsurf)%nodes)
167  ndlocal(1) = fstrsolid%contacts(i)%slave(j)
168  ndlocal(2:nnode+1) = fstrsolid%contacts(i)%master(ctsurf)%nodes(1:nnode)
169 
170  ! For CONTACTNEAR damping (especially S-Lagrange + frictionless),
171  ! we still need slave-master connectivity to assemble damping terms.
172  necessary_to_insert_node_pair = necessary_to_insert_node .or. is_damping_active_flag
173 
174  if( is_contact_active_flag ) then
175  do k=1,nlag
176  if( contact_algo == kcaslagrange ) count_lagrange = count_lagrange + 1
177  call hecmw_ass_noderelated_from_contact_pair(np, nnode, ndlocal, count_lagrange, permission, &
178  & necessary_to_insert_node_pair, list_noderelated_org, list_noderelated, countnon0lu_node, countnon0lu_lagrange )
179  enddo
180  else
181  ! NEAR damping only: no Lagrange multiplier, insert connectivity once
182  call hecmw_ass_noderelated_from_contact_pair(np, nnode, ndlocal, 0, permission, &
183  & necessary_to_insert_node_pair, list_noderelated_org, list_noderelated, countnon0lu_node, countnon0lu_lagrange )
184  endif
185 
186  end if
187 
188  enddo
189 
190  enddo
191 
192  do i = 1, fstrsolid%n_embeds
193 
194  grpid = fstrsolid%embeds(i)%group
195  if( .not. fstr_isembedactive( fstrsolid, grpid, cstep ) ) cycle
196 
197  necessary_to_insert_node = ( contact_algo == kcaalagrange )
198 
199  nlag = 3
200  if( contact_algo == kcaalagrange ) nlag = 1
201  permission = .true.
202 
203  do j = 1, size(fstrsolid%embeds(i)%slave)
204 
205  if( .not. is_contact_active(fstrsolid%embeds(i)%states(j)%state) ) cycle
206  ctsurf = fstrsolid%embeds(i)%states(j)%surface
207  etype = fstrsolid%embeds(i)%master(ctsurf)%etype
208  nnode = size(fstrsolid%embeds(i)%master(ctsurf)%nodes)
209  ndlocal(1) = fstrsolid%embeds(i)%slave(j)
210  ndlocal(2:nnode+1) = fstrsolid%embeds(i)%master(ctsurf)%nodes(1:nnode)
211 
212  do k=1,nlag
213  if( contact_algo == kcaslagrange ) count_lagrange = count_lagrange + 1
214  call hecmw_ass_noderelated_from_contact_pair(np, nnode, ndlocal, count_lagrange, permission, &
215  & necessary_to_insert_node, list_noderelated_org, list_noderelated, countnon0lu_node, countnon0lu_lagrange )
216  enddo
217  enddo
218 
219  enddo
220 
221  end subroutine getnewlistofrelatednodesandlagrangemultipliers
222 
224  subroutine fstr_copy_lagrange_contact(fstrSOLID,hecLagMAT)
225 
226  type(fstr_solid) :: fstrsolid
227  type(hecmwst_matrix_lagrange) :: heclagmat
228  integer (kind=kint) :: id_lagrange, algtype, i, j, k, nlag, slave_node
229 
230  id_lagrange = 0
231 
232  do i = 1, fstrsolid%n_contacts
233 
234  algtype = fstrsolid%contacts(i)%algtype
235  nlag = fstr_get_num_lagrange_pernode(algtype)
236 
237  do j = 1, size(fstrsolid%contacts(i)%slave)
238  if( .not. is_contact_active(fstrsolid%contacts(i)%states(j)%state) ) cycle
239  slave_node = fstrsolid%contacts(i)%slave(j)
240  heclagmat%lag_node_table(slave_node) = id_lagrange + 1
241  do k=1,nlag
242  id_lagrange = id_lagrange + 1
243  heclagmat%Lagrange(id_lagrange)=fstrsolid%contacts(i)%states(j)%multiplier(k)
244  enddo
245  enddo
246  enddo
247 
248  do i = 1, fstrsolid%n_embeds
249  nlag = 3
250  do j = 1, size(fstrsolid%embeds(i)%slave)
251  if( .not. is_contact_active(fstrsolid%embeds(i)%states(j)%state) ) cycle
252  slave_node = fstrsolid%embeds(i)%slave(j)
253  heclagmat%lag_node_table(slave_node) = id_lagrange + 1
254  do k=1,nlag
255  id_lagrange = id_lagrange + 1
256  heclagmat%Lagrange(id_lagrange)=fstrsolid%embeds(i)%states(j)%multiplier(k)
257  enddo
258  enddo
259  enddo
260 
261  end subroutine fstr_copy_lagrange_contact
262 
264  logical function fstr_is_matrixstruct_symmetric(fstrSOLID,hecMESH)
265 
266  type(fstr_solid ) :: fstrsolid
267  type(hecmwst_local_mesh) :: hecmesh
268  integer (kind=kint) :: is_in_contact
269 
270  is_in_contact = 0
271  if( fstrsolid%n_contacts>0 ) then
272  if( any(fstrsolid%contacts(:)%fcoeff /= 0.0d0) ) is_in_contact = 1
273  endif
274  call hecmw_allreduce_i1(hecmesh, is_in_contact, hecmw_max)
275  if( is_in_contact == 0 .and. hecmesh%n_dof /= 4 ) then
277  else
279  endif
280 
281  end function fstr_is_matrixstruct_symmetric
282 
283 
284 end module fstr_matrix_con_contact
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
integer, parameter fe_tri3n
Definition: element.f90:69
integer, parameter fe_quad4n
Definition: element.f90:72
This module provides functions of reconstructing.
subroutine, public fstr_mat_con_contact(cstep, contact_algo, hecMAT, fstrSOLID, hecLagMAT, infoCTChange, conMAT, is_contact_active_flag)
this subroutine reconstructs node-based (stiffness) matrix structure \corresponding to contact state
integer(kind=kint) function, public fstr_get_num_lagrange_pernode(algtype)
subroutine, public fstr_save_originalmatrixstructure(hecMAT)
This subroutine saves original matrix structure constructed originally by hecMW_matrix.
logical function, public fstr_is_matrixstruct_symmetric(fstrSOLID, hecMESH)
this function judges whether sitiffness matrix is symmetric or not
Contact damping module for CONTACTNEAR state.
pure logical function, public is_damping_enabled(contact)
Check if damping is enabled for a contact pair.
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
integer(kind=kint), parameter kcaalagrange
Definition: m_fstr.F90:62