FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
surf_ele.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 !-------------------------------------------------------------------------------
9  use hecmw_util
10  use elementinfo
11  use m_utilities
12  use bucket_search
13 
14  implicit none
15 
16  integer(kind=kint), parameter :: l_max_surface_node =20
17  integer(kind=kint), parameter :: l_max_elem_node = 100
18  integer(kind=kint), parameter :: l_max_elem_surf = 6
19 
20  integer(kind=kint), parameter :: n_neighbor_max_init = 8
21 
24  integer(kind=kint) :: eid
25  integer(kind=kint) :: etype
26  integer(kind=kint), pointer :: nodes(:)=>null()
27  integer(kind=kint) :: n_neighbor
28  integer(kind=kint) :: n_neighbor_max
29  integer(kind=kint), pointer :: neighbor(:)=>null()
30  real(kind=kreal) :: reflen
31  real(kind=kreal) :: xavg(3)
32  real(kind=kreal) :: dmax
33  integer(kind=kint) :: bktid
34  end type tsurfelement
35 
36  integer(kind=kint), parameter, private :: debug = 0
37 
38 contains
39 
41  subroutine initialize_surf( eid, etype, nsurf, surf )
42  integer(kind=kint), intent(in) :: eid
43  integer(kind=kint), intent(in) :: etype
44  integer(kind=kint), intent(in) :: nsurf
45  type(tsurfelement), intent(inout) :: surf
46  integer(kind=kint) :: i, n, outtype, nodes(100)
47  surf%eid = eid
48 
49  if( nsurf > 0 ) then
50  call getsubface( etype, nsurf, outtype, nodes )
51  surf%etype = outtype
52  n=getnumberofnodes( outtype )
53  else !treat 3d master element as it is
54  surf%etype = etype
55  n=getnumberofnodes( etype )
56  do i=1,n
57  nodes(i) = i
58  enddo
59  endif
60  allocate( surf%nodes(n) )
61  surf%nodes(1:n)=nodes(1:n)
62  surf%n_neighbor = 0
63  surf%n_neighbor_max = 0
64  surf%reflen = -1.d0
65  surf%xavg(:) = 0.d0
66  surf%dmax = -1.d0
67  surf%bktID = -1
68  end subroutine
69 
71  subroutine finalize_surf( surf )
72  type(tsurfelement), intent(inout) :: surf
73  if( associated(surf%nodes) ) deallocate( surf%nodes )
74  if( associated(surf%neighbor) ) deallocate( surf%neighbor )
75  end subroutine
76 
78  subroutine write_surf( file, surf )
79  integer(kind=kint), intent(in) :: file
80  type(tsurfelement), intent(in) :: surf
81  integer(kind=kint) :: i
82  write(file,*) "Element:",surf%eid,"Surface type=",surf%etype
83  if( associated( surf%nodes ) ) &
84  write(file,*) ( surf%nodes(i),i=1,size(surf%nodes) )
85  if( associated( surf%neighbor ) ) &
86  write(file,*) ( surf%neighbor(i),i=1,surf%n_neighbor )
87  end subroutine
88 
90  subroutine find_surface_neighbor( surf, bktDB )
91  type(tsurfelement), intent(inout) :: surf(:)
92  type(bucketdb), intent(in) :: bktDB
93  integer(kind=kint) :: i, j, ii,jj, nd1, nd2, nsurf
94  integer(kind=kint) :: k, oldsize, newsize
95  integer(kind=kint), pointer :: dumarray(:) => null()
96  integer(kind=kint) :: bktID, ncand, js
97  integer(kind=kint), allocatable :: indexSurf(:)
98  if (debug >= 1) write(0,*) 'DEBUG: find_surface_neighbor: start'
99 
100  nsurf = size(surf)
101 
102  !$omp parallel do default(none), &
103  !$omp&private(i,ii,nd1,j,jj,nd2,oldsize,newsize,dumarray,k,bktID,ncand,indexSurf,js), &
104  !$omp&shared(nsurf,surf,bktDB)
105  do i=1,nsurf
106  bktid = bucketdb_getbucketid(bktdb, surf(i)%xavg)
107  ncand = bucketdb_getnumcand(bktdb, bktid)
108  if (ncand == 0) cycle
109  allocate(indexsurf(ncand))
110  call bucketdb_getcand(bktdb, bktid, ncand, indexsurf)
111  jloop: do js=1,ncand
112  j = indexsurf(js)
113  if( i==j ) cycle
114  if( associated(surf(i)%neighbor) ) then
115  if ( any( surf(i)%neighbor(1:surf(i)%n_neighbor)==j ) ) cycle
116  endif
117  do ii=1, size(surf(i)%nodes)
118  nd1 = surf(i)%nodes(ii)
119  do jj=1, size(surf(j)%nodes)
120  nd2 = surf(j)%nodes(jj)
121  if( nd1==nd2 ) then
122  surf(i)%n_neighbor = surf(i)%n_neighbor+1
123 
124  if( .not. associated(surf(i)%neighbor) ) then
125  allocate( surf(i)%neighbor(n_neighbor_max_init) )
126  surf(i)%n_neighbor_max = n_neighbor_max_init
127  surf(i)%neighbor(1) = j
128  else if( surf(i)%n_neighbor > surf(i)%n_neighbor_max ) then
129  oldsize = surf(i)%n_neighbor_max
130  newsize = oldsize * 2
131  dumarray => surf(i)%neighbor
132  allocate( surf(i)%neighbor(newsize) )
133  surf(i)%n_neighbor_max = newsize
134  do k=1,oldsize
135  surf(i)%neighbor(k) = dumarray(k)
136  enddo
137  surf(i)%neighbor(oldsize+1) = j
138  deallocate( dumarray )
139  else
140  surf(i)%neighbor(surf(i)%n_neighbor) = j
141  endif
142 
143  cycle jloop
144  endif
145  enddo
146  enddo
147  enddo jloop
148  deallocate(indexsurf)
149  enddo
150  !$omp end parallel do
151 
152  if (debug >= 1) write(0,*) 'DEBUG: find_surface_neighbor: end'
153  end subroutine
154 
156  integer(kind=kint) function next_position( surf, cpos )
157  type(tsurfelement), intent(in) :: surf
158  real(kind=kreal), intent(in) :: cpos(2)
159 
160  integer(kind=kint) :: i
161  real(kind=kreal) :: maxv(3)
162  next_position = surf%eid
163  if( .not. associated(surf%neighbor) ) return ! do nothing when not initialized
164  maxv(:) = 0.d0
165  select case(surf%etype)
166  case( fe_tri3n, fe_tri6n )
167  if( all(cpos>0.d0) .and. all(cpos<1.d0) ) return
168  if( size(surf%neighbor)<3 ) return
169  if( all(cpos(:)>0.d0) ) return
170  do i=1,3
171  if( cpos(i)< 0.d0 ) maxv(i) = -cpos(i)
172  enddo
173  next_position = maxloc(maxv,1)
174  case( fe_quad4n, fe_quad8n )
175  if( all(cpos>-1.d0) .and. all(cpos<1.d0) ) return
176  if( size(surf%neighbor)<4 ) return
177  if( cpos(1)<-1.d0 .and. dabs(cpos(2))<1.d0 ) then
178  next_position = 1
179  elseif( cpos(1)>1.d0 .and. dabs(cpos(2))<1.d0 ) then
180  next_position = 3
181  elseif( dabs(cpos(1))<1.d0 .and. cpos(2)<-1.d0 ) then
182  next_position = 2
183  elseif( dabs(cpos(1))<1.d0 .and. cpos(2)>1.d0 ) then
184  next_position = 4
185  elseif( cpos(1)<-1.d0 .and. cpos(2)<-1.d0 ) then
186  if( cpos(1)>cpos(2) ) then
187  next_position = 2
188  else
189  next_position = 1
190  endif
191  elseif( cpos(1)<-1.d0 .and. cpos(2)>1.d0 ) then
192  if( dabs(cpos(1))>cpos(2) ) then
193  next_position = 1
194  else
195  next_position = 4
196  endif
197  elseif( cpos(1)>1.d0 .and. cpos(2)<-1.d0 ) then
198  if( cpos(1)>dabs(cpos(2)) ) then
199  next_position = 3
200  else
201  next_position = 2
202  endif
203  elseif( cpos(1)>1.d0 .and. cpos(2)>1.d0 ) then
204  if( cpos(1)>cpos(2) ) then
205  next_position = 3
206  else
207  next_position = 4
208  endif
209  endif
210  case default
211  stop "type of surface element not defined"
212  end select
213  next_position = surf%neighbor( next_position )
214 
215  end function next_position
216 
218  subroutine update_surface_reflen( surf, coord )
219  type(tsurfelement), intent(inout) :: surf(:)
220  real(kind=kreal), intent(in) :: coord(:)
221  real(kind=kreal) :: elem(3, l_max_surface_node), r0(2)
222  integer(kind=kint) :: nn, i, j, iss
223  !$omp parallel do default(none) private(j,nn,i,iss,elem,r0) shared(surf,coord)
224  do j=1,size(surf)
225  nn = size(surf(j)%nodes)
226  do i=1,nn
227  iss = surf(j)%nodes(i)
228  elem(1:3,i) = coord(3*iss-2:3*iss)
229  enddo
230  call getelementcenter( surf(j)%etype, r0 )
231  surf(j)%reflen = getreferencelength( surf(j)%etype, nn, r0, elem )
232  enddo
233  !$omp end parallel do
234  end subroutine update_surface_reflen
235 
237  subroutine update_surface_box_info( surf, currpos )
238  type(tsurfelement), intent(inout) :: surf(:)
239  real(kind=kreal), intent(in) :: currpos(:)
240  integer(kind=kint) :: nn, i, j, iss
241  real(kind=kreal) :: elem(3,l_max_surface_node),xmin(3), xmax(3), xsum(3), lx(3)
242  !$omp parallel do default(none) private(i,nn,iss,elem,xmin,xmax,xsum,lx) shared(surf,currpos)
243  do j=1, size(surf)
244  nn = size(surf(j)%nodes)
245  do i=1,nn
246  iss = surf(j)%nodes(i)
247  elem(1:3,i) = currpos(3*iss-2:3*iss)
248  enddo
249  do i=1,3
250  xmin(i) = minval(elem(i,1:nn))
251  xmax(i) = maxval(elem(i,1:nn))
252  xsum(i) = sum(elem(i,1:nn))
253  enddo
254  surf(j)%xavg(:) = xsum(:) / nn
255  lx(:) = xmax(:) - xmin(:)
256  surf(j)%dmax = maxval(lx) * 0.5d0
257  enddo
258  !$omp end parallel do
259  end subroutine update_surface_box_info
260 
262  logical function is_in_surface_box(surf, x0, exp_rate)
263  type(tsurfelement), intent(inout) :: surf
264  real(kind=kreal), intent(in) :: x0(3)
265  real(kind=kreal), intent(in) :: exp_rate
266  real(kind=kreal) :: dist(3), er
267  er = max(exp_rate, 1.d0)
268  dist(:) = abs(x0(:) - surf%xavg(:))
269  if ( maxval(dist(:)) < surf%dmax * er ) then
270  is_in_surface_box = .true.
271  else
272  is_in_surface_box = .false.
273  endif
274  end function is_in_surface_box
275 
277  subroutine update_surface_bucket_info(surf, bktDB)
278  type(tsurfelement), intent(inout) :: surf(:)
279  type(bucketdb), intent(inout) :: bktDB
280  real(kind=kreal) :: x_min(3), x_max(3), d_max
281  integer(kind=kint) :: nsurf, i, j
282  if (debug >= 1) write(0,*) 'DEBUG: update_surface_bucket_info: start'
283  nsurf = size(surf)
284  if (nsurf == 0) return
285  x_min(:) = surf(1)%xavg(:)
286  x_max(:) = surf(1)%xavg(:)
287  d_max = surf(1)%dmax*2.d0
288  do i = 2, nsurf
289  do j = 1, 3
290  if (surf(i)%xavg(j) < x_min(j)) x_min(j) = surf(i)%xavg(j)
291  if (surf(i)%xavg(j) > x_max(j)) x_max(j) = surf(i)%xavg(j)
292  if (surf(i)%dmax > d_max) d_max = surf(i)%dmax
293  enddo
294  enddo
295  do j = 1, 3
296  x_min(j) = x_min(j) - d_max
297  x_max(j) = x_max(j) + d_max
298  enddo
299  call bucketdb_setup(bktdb, x_min, x_max, d_max, nsurf)
300  !$omp parallel do default(none) private(i) shared(nsurf,surf,bktDB)
301  do i = 1, nsurf
302  surf(i)%bktID = bucketdb_getbucketid(bktdb, surf(i)%xavg)
303  call bucketdb_registerpre(bktdb, surf(i)%bktID)
304  enddo
305  !$omp end parallel do
306  call bucketdb_allocate(bktdb)
307  do i = 1, nsurf
308  call bucketdb_register(bktdb, surf(i)%bktID, i)
309  enddo
310  if (debug >= 1) write(0,*) 'DEBUG: update_surface_bucket_info: end'
311  end subroutine update_surface_bucket_info
312 
313 end module msurfelement
This module provides bucket-search functionality It provides definition of bucket info and its access...
subroutine, public bucketdb_allocate(bktdb)
Allocate memory before actually registering members Before allocating memory, bucketDB_registerPre ha...
subroutine, public bucketdb_registerpre(bktdb, bid)
Pre-register for just counting members to be actually registered Bucket ID has to be obtained with bu...
integer(kind=kint) function, public bucketdb_getbucketid(bktdb, x)
Get bucket ID that includes given point.
subroutine, public bucketdb_register(bktdb, bid, sid)
Register member Before actually register, bucketDB_allocate has to be called.
integer(kind=kint) function, public bucketdb_getnumcand(bktdb, bid)
Get number of candidates within neighboring buckets of a given bucket Bucket ID has to be obtained wi...
subroutine, public bucketdb_getcand(bktdb, bid, ncand, cand)
Get candidates within neighboring buckets of a given bucket Number of candidates has to be obtained w...
subroutine, public bucketdb_setup(bktdb, x_min, x_max, dmin, n_tot)
Setup basic info of buckets.
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
integer(kind=kind(2)) function getnumberofnodes(etype)
Obtain number of nodes of the element.
Definition: element.f90:131
subroutine getelementcenter(fetype, localcoord)
Return natural coordinate of the center of element.
Definition: element.f90:1012
integer, parameter fe_tri6n
Definition: element.f90:70
subroutine getsubface(intype, innumber, outtype, nodes)
Find the definition of surface of the element.
Definition: element.f90:193
integer, parameter fe_tri3n
Definition: element.f90:69
integer, parameter fe_quad4n
Definition: element.f90:72
real(kind=kreal) function getreferencelength(fetype, nn, localcoord, elecoord)
This function calculates reference length at a point in surface.
Definition: element.f90:1266
integer, parameter fe_quad8n
Definition: element.f90:73
I/O and Utility.
Definition: hecmw_util_f.F90:7
This module provides aux functions.
Definition: utilities.f90:6
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
integer(kind=kint), parameter l_max_elem_node
Definition: surf_ele.f90:17
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
integer(kind=kint), parameter l_max_elem_surf
Definition: surf_ele.f90:18
subroutine find_surface_neighbor(surf, bktDB)
Find neighboring surface elements.
Definition: surf_ele.f90:91
integer(kind=kint), parameter n_neighbor_max_init
Definition: surf_ele.f90:20
subroutine update_surface_bucket_info(surf, bktDB)
Update bucket info for searching surface elements.
Definition: surf_ele.f90:278
integer(kind=kint), parameter l_max_surface_node
Definition: surf_ele.f90:16
subroutine finalize_surf(surf)
Memory management subroutine.
Definition: surf_ele.f90:72
logical function is_in_surface_box(surf, x0, exp_rate)
Check if given point is within cubic box including surface element.
Definition: surf_ele.f90:263
integer(kind=kint) function next_position(surf, cpos)
Tracing next contact position.
Definition: surf_ele.f90:157
Structure to define surface group.
Definition: surf_ele.f90:23