FrontISTR  5.7.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, only: kint, kreal
10  implicit none
11 
12  integer(kind=kint), parameter :: l_max_surface_node =20
13  integer(kind=kint), parameter :: l_max_elem_node = 100
14  integer(kind=kint), parameter :: l_max_elem_surf = 6
15 
16  integer(kind=kint), parameter :: n_neighbor_max_init = 8
17 
20  integer(kind=kint) :: eid
21  integer(kind=kint) :: etype
22  integer(kind=kint), pointer :: nodes(:)=>null()
23  integer(kind=kint) :: n_neighbor
24  integer(kind=kint) :: n_neighbor_max
25  integer(kind=kint), pointer :: neighbor(:)=>null()
26  real(kind=kreal) :: reflen
27  real(kind=kreal) :: xavg(3)
28  real(kind=kreal) :: dmax
29  integer(kind=kint) :: bktid
30  end type tsurfelement
31 
32  integer(kind=kint), parameter, private :: debug = 0
33 
34 contains
35 
37  subroutine initialize_surf( eid, etype, nsurf, surf )
39  integer(kind=kint), intent(in) :: eid
40  integer(kind=kint), intent(in) :: etype
41  integer(kind=kint), intent(in) :: nsurf
42  type(tsurfelement), intent(inout) :: surf
43  integer(kind=kint) :: i, n, outtype, nodes(100)
44  surf%eid = eid
45 
46  if( nsurf > 0 ) then
47  call getsubface( etype, nsurf, outtype, nodes )
48  surf%etype = outtype
49  n=getnumberofnodes( outtype )
50  else !treat 3d master element as it is
51  surf%etype = etype
52  n=getnumberofnodes( etype )
53  do i=1,n
54  nodes(i) = i
55  enddo
56  endif
57  allocate( surf%nodes(n) )
58  surf%nodes(1:n)=nodes(1:n)
59  surf%n_neighbor = 0
60  surf%n_neighbor_max = 0
61  surf%reflen = -1.d0
62  surf%xavg(:) = 0.d0
63  surf%dmax = -1.d0
64  surf%bktID = -1
65  end subroutine
66 
68  subroutine finalize_surf( surf )
69  type(tsurfelement), intent(inout) :: surf
70  if( associated(surf%nodes) ) deallocate( surf%nodes )
71  if( associated(surf%neighbor) ) deallocate( surf%neighbor )
72  end subroutine
73 
75  subroutine write_surf( file, surf )
76  integer(kind=kint), intent(in) :: file
77  type(tsurfelement), intent(in) :: surf
78  integer(kind=kint) :: i
79  write(file,*) "Element:",surf%eid,"Surface type=",surf%etype
80  if( associated( surf%nodes ) ) &
81  write(file,*) ( surf%nodes(i),i=1,size(surf%nodes) )
82  if( associated( surf%neighbor ) ) &
83  write(file,*) ( surf%neighbor(i),i=1,surf%n_neighbor )
84  end subroutine
85 
87  subroutine find_surface_neighbor( surf, bktDB )
89  use hecmw_util
90  use bucket_search
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 )
158  type(tsurfelement), intent(in) :: surf
159  real(kind=kreal), intent(in) :: cpos(2)
160 
161  integer(kind=kint) :: i
162  real(kind=kreal) :: maxv(3)
163  next_position = surf%eid
164  if( .not. associated(surf%neighbor) ) return ! do nothing when not initialized
165  maxv(:) = 0.d0
166  select case(surf%etype)
167  case( fe_tri3n, fe_tri6n )
168  if( all(cpos>0.d0) .and. all(cpos<1.d0) ) return
169  if( size(surf%neighbor)<3 ) return
170  if( all(cpos(:)>0.d0) ) return
171  do i=1,3
172  if( cpos(i)< 0.d0 ) maxv(i) = -cpos(i)
173  enddo
174  next_position = maxloc(maxv,1)
175  case( fe_quad4n, fe_quad8n )
176  if( all(cpos>-1.d0) .and. all(cpos<1.d0) ) return
177  if( size(surf%neighbor)<4 ) return
178  if( cpos(1)<-1.d0 .and. dabs(cpos(2))<1.d0 ) then
179  next_position = 1
180  elseif( cpos(1)>1.d0 .and. dabs(cpos(2))<1.d0 ) then
181  next_position = 3
182  elseif( dabs(cpos(1))<1.d0 .and. cpos(2)<-1.d0 ) then
183  next_position = 2
184  elseif( dabs(cpos(1))<1.d0 .and. cpos(2)>1.d0 ) then
185  next_position = 4
186  elseif( cpos(1)<-1.d0 .and. cpos(2)<-1.d0 ) then
187  if( cpos(1)>cpos(2) ) then
188  next_position = 2
189  else
190  next_position = 1
191  endif
192  elseif( cpos(1)<-1.d0 .and. cpos(2)>1.d0 ) then
193  if( dabs(cpos(1))>cpos(2) ) then
194  next_position = 1
195  else
196  next_position = 4
197  endif
198  elseif( cpos(1)>1.d0 .and. cpos(2)<-1.d0 ) then
199  if( cpos(1)>dabs(cpos(2)) ) then
200  next_position = 3
201  else
202  next_position = 2
203  endif
204  elseif( cpos(1)>1.d0 .and. cpos(2)>1.d0 ) then
205  if( cpos(1)>cpos(2) ) then
206  next_position = 3
207  else
208  next_position = 4
209  endif
210  endif
211  case default
212  stop "type of surface element not defined"
213  end select
214  next_position = surf%neighbor( next_position )
215 
216  end function next_position
217 
219  subroutine update_surface_reflen( surf, coord )
221  type(tsurfelement), intent(inout) :: surf(:)
222  real(kind=kreal), intent(in) :: coord(:)
223  real(kind=kreal) :: elem(3, l_max_surface_node), r0(2)
224  integer(kind=kint) :: nn, i, j, iss
225  !$omp parallel do default(none) private(j,nn,i,iss,elem,r0) shared(surf,coord)
226  do j=1,size(surf)
227  nn = size(surf(j)%nodes)
228  do i=1,nn
229  iss = surf(j)%nodes(i)
230  elem(1:3,i) = coord(3*iss-2:3*iss)
231  enddo
232  call getelementcenter( surf(j)%etype, r0 )
233  surf(j)%reflen = getreferencelength( surf(j)%etype, nn, r0, elem )
234  enddo
235  !$omp end parallel do
236  end subroutine update_surface_reflen
237 
239  subroutine update_surface_box_info( surf, currpos )
240  type(tsurfelement), intent(inout) :: surf(:)
241  real(kind=kreal), intent(in) :: currpos(:)
242  integer(kind=kint) :: nn, i, j, iss
243  real(kind=kreal) :: elem(3,l_max_surface_node),xmin(3), xmax(3), xsum(3), lx(3)
244  !$omp parallel do default(none) private(i,nn,iss,elem,xmin,xmax,xsum,lx) shared(surf,currpos)
245  do j=1, size(surf)
246  nn = size(surf(j)%nodes)
247  do i=1,nn
248  iss = surf(j)%nodes(i)
249  elem(1:3,i) = currpos(3*iss-2:3*iss)
250  enddo
251  do i=1,3
252  xmin(i) = minval(elem(i,1:nn))
253  xmax(i) = maxval(elem(i,1:nn))
254  xsum(i) = sum(elem(i,1:nn))
255  enddo
256  surf(j)%xavg(:) = xsum(:) / nn
257  lx(:) = xmax(:) - xmin(:)
258  surf(j)%dmax = maxval(lx) * 0.5d0
259  enddo
260  !$omp end parallel do
261  end subroutine update_surface_box_info
262 
264  logical function is_in_surface_box(surf, x0, exp_rate)
265  type(tsurfelement), intent(inout) :: surf
266  real(kind=kreal), intent(in) :: x0(3)
267  real(kind=kreal), intent(in) :: exp_rate
268  real(kind=kreal) :: dist(3), er
269  er = max(exp_rate, 1.d0)
270  dist(:) = abs(x0(:) - surf%xavg(:))
271  if ( maxval(dist(:)) < surf%dmax * er ) then
272  is_in_surface_box = .true.
273  else
274  is_in_surface_box = .false.
275  endif
276  end function is_in_surface_box
277 
279  subroutine update_surface_bucket_info(surf, bktDB)
281  type(tsurfelement), intent(inout) :: surf(:)
282  type(bucketdb), intent(inout) :: bktDB
283  real(kind=kreal) :: x_min(3), x_max(3), d_max
284  integer(kind=kint) :: nsurf, i, j
285  if (debug >= 1) write(0,*) 'DEBUG: update_surface_bucket_info: start'
286  nsurf = size(surf)
287  if (nsurf == 0) return
288  x_min(:) = surf(1)%xavg(:)
289  x_max(:) = surf(1)%xavg(:)
290  d_max = surf(1)%dmax*2.d0
291  do i = 2, nsurf
292  do j = 1, 3
293  if (surf(i)%xavg(j) < x_min(j)) x_min(j) = surf(i)%xavg(j)
294  if (surf(i)%xavg(j) > x_max(j)) x_max(j) = surf(i)%xavg(j)
295  if (surf(i)%dmax > d_max) d_max = surf(i)%dmax
296  enddo
297  enddo
298  do j = 1, 3
299  x_min(j) = x_min(j) - d_max
300  x_max(j) = x_max(j) + d_max
301  enddo
302  call bucketdb_setup(bktdb, x_min, x_max, d_max, nsurf)
303  !$omp parallel do default(none) private(i) shared(nsurf,surf,bktDB)
304  do i = 1, nsurf
305  surf(i)%bktID = bucketdb_getbucketid(bktdb, surf(i)%xavg)
306  call bucketdb_registerpre(bktdb, surf(i)%bktID)
307  enddo
308  !$omp end parallel do
309  call bucketdb_allocate(bktdb)
310  do i = 1, nsurf
311  call bucketdb_register(bktdb, surf(i)%bktID, i)
312  enddo
313  if (debug >= 1) write(0,*) 'DEBUG: update_surface_bucket_info: end'
314  end subroutine update_surface_bucket_info
315 
316 end module msurfelement
msurfelement::n_neighbor_max_init
integer(kind=kint), parameter n_neighbor_max_init
Definition: surf_ele.f90:16
elementinfo::fe_quad4n
integer, parameter fe_quad4n
Definition: element.f90:72
msurfelement::l_max_elem_surf
integer(kind=kint), parameter l_max_elem_surf
Definition: surf_ele.f90:14
m_utilities
This module provides aux functions.
Definition: utilities.f90:6
msurfelement::next_position
integer(kind=kint) function next_position(surf, cpos)
Tracing next contact position.
Definition: surf_ele.f90:157
msurfelement::finalize_surf
subroutine finalize_surf(surf)
Memory management subroutine.
Definition: surf_ele.f90:69
bucket_search::bucketdb_registerpre
subroutine, public bucketdb_registerpre(bktdb, bid)
Pre-register for just counting members to be actually registered Bucket ID has to be obtained with bu...
Definition: bucket_search.f90:265
elementinfo::fe_tri3n
integer, parameter fe_tri3n
Definition: element.f90:69
elementinfo::fe_quad8n
integer, parameter fe_quad8n
Definition: element.f90:73
bucket_search
This module provides bucket-search functionality It provides definition of bucket info and its access...
Definition: bucket_search.f90:7
msurfelement::update_surface_bucket_info
subroutine update_surface_bucket_info(surf, bktDB)
Update bucket info for searching surface elements.
Definition: surf_ele.f90:280
msurfelement::write_surf
subroutine write_surf(file, surf)
Write out elemental surface.
Definition: surf_ele.f90:76
bucket_search::bucketdb_register
subroutine, public bucketdb_register(bktdb, bid, sid)
Register member Before actually register, bucketDB_allocate has to be called.
Definition: bucket_search.f90:299
msurfelement::is_in_surface_box
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:265
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
msurfelement::tsurfelement
Structure to define surface group.
Definition: surf_ele.f90:19
msurfelement::update_surface_box_info
subroutine update_surface_box_info(surf, currpos)
Update info of cubic box including surface elements.
Definition: surf_ele.f90:240
msurfelement::find_surface_neighbor
subroutine find_surface_neighbor(surf, bktDB)
Find neighboring surface elements.
Definition: surf_ele.f90:88
elementinfo::getreferencelength
real(kind=kreal) function getreferencelength(fetype, nn, localcoord, elecoord)
This function calculates reference length at a point in surface.
Definition: element.f90:1259
msurfelement::l_max_elem_node
integer(kind=kint), parameter l_max_elem_node
Definition: surf_ele.f90:13
msurfelement::update_surface_reflen
subroutine update_surface_reflen(surf, coord)
Compute reference length of surface elements.
Definition: surf_ele.f90:220
elementinfo
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
bucket_search::bucketdb_getnumcand
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...
Definition: bucket_search.f90:313
bucket_search::bucketdb_getcand
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...
Definition: bucket_search.f90:344
hecmw
Definition: hecmw.f90:6
bucket_search::bucketdb_getbucketid
integer(kind=kint) function, public bucketdb_getbucketid(bktdb, x)
Get bucket ID that includes given point.
Definition: bucket_search.f90:244
elementinfo::getsubface
subroutine getsubface(intype, innumber, outtype, nodes)
Find the definition of surface of the element.
Definition: element.f90:193
bucket_search::bucketdb_setup
subroutine, public bucketdb_setup(bktdb, x_min, x_max, dmin, n_tot)
Setup basic info of buckets.
Definition: bucket_search.f90:175
msurfelement::l_max_surface_node
integer(kind=kint), parameter l_max_surface_node
Definition: surf_ele.f90:12
msurfelement::initialize_surf
subroutine initialize_surf(eid, etype, nsurf, surf)
Initializer.
Definition: surf_ele.f90:38
elementinfo::fe_tri6n
integer, parameter fe_tri6n
Definition: element.f90:70
elementinfo::getelementcenter
subroutine getelementcenter(fetype, localcoord)
Return natural coordinate of the center of surface element.
Definition: element.f90:1012
elementinfo::getnumberofnodes
integer(kind=kind(2)) function getnumberofnodes(etype)
Obtain number of nodes of the element.
Definition: element.f90:131
bucket_search::bucketdb_allocate
subroutine, public bucketdb_allocate(bktdb)
Allocate memory before actually registering members Before allocating memory, bucketDB_registerPre ha...
Definition: bucket_search.f90:278
msurfelement
This module manages surface elements in 3D It provides basic definition of surface elements (triangla...
Definition: surf_ele.f90:8