18 public :: hecmwst_matrix_lagrange
23 integer(kind=kint),
save :: NPL_org, NPU_org
24 type(nodeRelated),
pointer,
save :: list_nodeRelated_org(:) => null()
26 type(nodeRelated),
pointer :: list_nodeRelated(:) => null()
28 logical :: permission = .false.
33 integer(kind=kint) :: algtype
34 if( algtype == contactsslid .or. algtype == contactfslid )
then
36 else if( algtype == contacttied )
then
44 type(hecmwst_matrix) :: hecmat
46 if(
associated(list_noderelated_org) )
return
47 call hecmw_construct_noderelated_from_hecmat(hecmat, npl_org, npu_org, list_noderelated_org)
53 subroutine fstr_mat_con_contact(cstep,contact_algo,hecMAT,fstrSOLID,hecLagMAT,infoCTChange,conMAT,is_contact_active_flag)
55 integer(kind=kint) :: cstep
56 integer(kind=kint) :: contact_algo
57 type(hecmwst_matrix) :: hecmat
59 type(hecmwst_matrix_lagrange) :: heclagmat
60 type(fstr_info_contactchange) :: infoctchange
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
69 integer(kind=kint) :: i, j, grpid
70 integer(kind=kint) :: nlag
74 do i = 1, fstrsolid%n_contacts
75 grpid = fstrsolid%contacts(i)%group
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
84 do i = 1, fstrsolid%n_embeds
85 grpid = fstrsolid%embeds(i)%group
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
96 call hecmw_init_noderelated_from_org(hecmat%NP,num_lagrange,is_contact_active_flag,list_noderelated_org,list_noderelated)
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)
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)
116 if( is_contact_active_flag .and. contact_algo ==
kcaslagrange ) &
117 call fstr_copy_lagrange_contact(fstrsolid,heclagmat)
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
128 integer(kind=kint),
intent(inout) :: countnon0lu_node, countnon0lu_lagrange
129 type(noderelated),
pointer,
intent(inout) :: list_noderelated(:)
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
140 do i = 1, fstrsolid%n_contacts
142 grpid = fstrsolid%contacts(i)%group
145 fcoeff = fstrsolid%contacts(i)%fcoeff
146 necessary_to_insert_node = ( fcoeff /= 0.0d0 .or. contact_algo ==
kcaalagrange )
148 algtype = fstrsolid%contacts(i)%algtype
151 if( algtype == contacttied ) permission = .true.
153 do j = 1,
size(fstrsolid%contacts(i)%slave)
155 is_contact_active_flag = is_contact_active(fstrsolid%contacts(i)%states(j)%state)
157 is_damping_active_flag = fstrsolid%contacts(i)%states(j)%state == contactnear .and. &
160 if( is_contact_active_flag .or. is_damping_active_flag )
then
162 ctsurf = fstrsolid%contacts(i)%states(j)%surface
163 etype = fstrsolid%contacts(i)%master(ctsurf)%etype
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)
172 necessary_to_insert_node_pair = necessary_to_insert_node .or. is_damping_active_flag
174 if( is_contact_active_flag )
then
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 )
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 )
192 do i = 1, fstrsolid%n_embeds
194 grpid = fstrsolid%embeds(i)%group
197 necessary_to_insert_node = ( contact_algo ==
kcaalagrange )
203 do j = 1,
size(fstrsolid%embeds(i)%slave)
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)
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 )
221 end subroutine getnewlistofrelatednodesandlagrangemultipliers
224 subroutine fstr_copy_lagrange_contact(fstrSOLID,hecLagMAT)
227 type(hecmwst_matrix_lagrange) :: heclagmat
228 integer (kind=kint) :: id_lagrange, algtype, i, j, k, nlag, slave_node
232 do i = 1, fstrsolid%n_contacts
234 algtype = fstrsolid%contacts(i)%algtype
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
242 id_lagrange = id_lagrange + 1
243 heclagmat%Lagrange(id_lagrange)=fstrsolid%contacts(i)%states(j)%multiplier(k)
248 do i = 1, fstrsolid%n_embeds
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
255 id_lagrange = id_lagrange + 1
256 heclagmat%Lagrange(id_lagrange)=fstrsolid%embeds(i)%states(j)%multiplier(k)
261 end subroutine fstr_copy_lagrange_contact
267 type(hecmwst_local_mesh) :: hecmesh
268 integer (kind=kint) :: is_in_contact
271 if( fstrsolid%n_contacts>0 )
then
272 if( any(fstrsolid%contacts(:)%fcoeff /= 0.0d0) ) is_in_contact = 1
274 call hecmw_allreduce_i1(hecmesh, is_in_contact, hecmw_max)
275 if( is_in_contact == 0 .and. hecmesh%n_dof /= 4 )
then
This module encapsulate the basic functions of all elements provide by this software.
integer, parameter fe_tri3n
integer, parameter fe_quad4n
This module defines common data and basic structures for analysis.
logical function fstr_isembedactive(fstrSOLID, nbc, cstep)
logical function fstr_iscontactactive(fstrSOLID, nbc, cstep)
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
integer(kind=kint), parameter kcaalagrange