Neko  0.9.0
A portable framework for high-order spectral element flow simulations
parmetis.F90
Go to the documentation of this file.
1 ! Copyright (c) 2021, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 module parmetis
35  use comm
36  use point
37  use utils
38  use num_types
39  use mesh_field
40  use mesh
41  use, intrinsic :: iso_c_binding
42  implicit none
43  private
44 
45  integer, private, parameter :: metis_ok = 1, metis_error = -4
46 
48 
49  interface
50  integer (c_int) function parmetis_v3_partmeshkway &
51  (elmdist, eptr, eind, elmwgt, wgtflag, numflag, ncon, &
52  ncommonnodes, nparts, tpwgts, ubvec, options, edgecut, part) &
53  bind(c, name='ParMETIS_V3_PartMeshKway_wrapper')
54  use, intrinsic :: iso_c_binding
55  implicit none
56  ! idx_t variables
57  type(c_ptr), value :: elmdist, eptr, eind, elmwgt, wgtflag, &
58  numflag, ncon, ncommonnodes, nparts, options, edgecut, part
59  ! real_t variables
60  type(c_ptr), value :: tpwgts, ubvec
61  end function parmetis_v3_partmeshkway
62  end interface
63 
64  interface
65  integer (c_int) function parmetis_v3_partgeom &
66  (vtxdist, ndims, xyz, part) &
67  bind(c, name='ParMETIS_V3_PartGeom_wrapper')
68  use, intrinsic :: iso_c_binding
69  implicit none
70  ! idx_t variables
71  type(c_ptr), value :: vtxdist, ndims, part
72  ! real_t variables
73  type(c_ptr), value :: xyz
74  end function parmetis_v3_partgeom
75  end interface
76 
77 !
78 ! Define data types depending on ParMETIS' configuration
79 !
80 ! parmetis_real converts between NEKO and ParMETIS' real representation
81 !
82 ! parmetis_idx converts between NEKO and ParMETIS' integer representation
83 ! neko_idx converts between ParMETIS and NEKO's integer representation
84 !
85 
86 #ifdef HAVE_PARMETIS
87 #ifdef HAVE_PARMETIS_REAL64
88 #define M_REAL c_double
89 #define parmetis_real(i) real((i), i8)
90 #else
91 #define M_REAL c_float
92 #define parmetis_real(i) real((i), 4)
93 #endif
94 #ifdef HAVE_PARMETIS_INT64
95 #define M_INT c_int64_t
96 #define parmetis_idx(i) int((i), i8)
97 #define neko_idx(i) int((i), 4)
98 #else
99 #define M_INT c_int32_t
100 #define parmetis_idx(i) (i)
101 #define neko_idx(i) (i)
102 #endif
103 #endif
104 
105 contains
106 
107 #ifdef HAVE_PARMETIS
108 
110  subroutine parmetis_partmeshkway(msh, parts, weights, nprts)
111  type(mesh_t), intent(inout) :: msh
112  type(mesh_fld_t), intent(inout) :: parts
113  type(mesh_fld_t), intent(in), optional :: weights
114  integer, intent(in), optional :: nprts
115  integer(kind=M_INT), target :: wgtflag, numflag, ncon, ncommonnodes
116  integer(kind=M_INT), target :: nparts, options(3), edgecut, rcode
117  real(kind=m_real), allocatable, target, dimension(:) :: tpwgts, ubvec
118  integer(kind=M_INT), allocatable, target, dimension(:) :: &
119  elmdist, eptr, eind, elmwgt, part
120  integer :: i, j, k, ierr
121 
125  numflag = 0
126  ncon = 1
127  ncommonnodes = 2**(msh%gdim - 1)
128  options(1) = 1
129  options(2) = 1
130  options(3) = 15
131  wgtflag = 2
132 
133  if (present(nprts)) then
134  nparts = nprts
135  else
136  nparts = pe_size
137  end if
138 
139  allocate(elmdist(0:pe_size), eptr(0:msh%nelv))
140  allocate(eind(0:(msh%nelv * msh%npts)), part(msh%nelv))
141  allocate(elmwgt(msh%nelv), tpwgts(ncon * nparts), ubvec(ncon))
142 
143  call parmetis_dist(elmdist, msh%nelv)
144 
145  if (present(weights)) then
146  call parmetis_wgt(msh, elmwgt, tpwgts, ubvec, nparts, ncon, weights)
147  else
148  call parmetis_wgt(msh, elmwgt, tpwgts, ubvec, nparts, ncon)
149  end if
150 
151  eptr(0) = 0
152  do i = 1, msh%nelv
153  eptr(i) = parmetis_idx(eptr(i - 1) + msh%npts)
154  end do
155 
156  k = 0
157  do i = 1, msh%nelv
158  do j = 1, msh%npts
159  eind(k) = parmetis_idx(msh%elements(i)%e%pts(j)%p%id() - 1)
160  k = k + 1
161  end do
162  end do
163 
164  rcode = parmetis_v3_partmeshkway(c_loc(elmdist), c_loc(eptr), c_loc(eind), &
165  c_loc(elmwgt), c_loc(wgtflag), c_loc(numflag), c_loc(ncon), &
166  c_loc(ncommonnodes), c_loc(nparts), c_loc(tpwgts), c_loc(ubvec),&
167  c_loc(options), c_loc(edgecut), c_loc(part))
168 
169  if (rcode .eq. metis_ok) then
170  call parmetis_mark_parts(parts, msh, part)
171  else
172  call neko_error(rcode)
173  end if
174 
175  deallocate(elmdist, eptr, eind, part, elmwgt, tpwgts, ubvec)
176 
177  end subroutine parmetis_partmeshkway
178 
181  subroutine parmetis_partgeom(msh, parts)
182  type(mesh_t), intent(inout) :: msh
183  type(mesh_fld_t), intent(inout) :: parts
184  integer(kind=M_INT), target :: ndims
185  real(kind=m_real), allocatable, target, dimension(:) :: xyz
186  integer(kind=M_INT), allocatable, target, dimension(:) :: vtxdist, part
187  type(point_t) :: c
188  integer :: i, j, ierr, rcode
189 
190  ndims = msh%gdim
191 
192  allocate(part(msh%nelv), xyz(ndims * msh%nelv))
193  allocate(vtxdist(0:pe_size))
194 
195  call parmetis_dist(vtxdist, msh%nelv)
196 
197  i = 1
198  do j = 1, msh%nelv
199  c = msh%elements(j)%e%centroid()
200  xyz(i) = parmetis_real(c%x(1))
201  xyz(i + 1) = parmetis_real(c%x(2))
202  xyz(i + 2) = parmetis_real(c%x(3))
203  i = i + 3
204  end do
205 
206  rcode = parmetis_v3_partgeom(c_loc(vtxdist), c_loc(ndims), &
207  c_loc(xyz), c_loc(part))
208 
209  if (rcode .eq. metis_ok) then
210  call parmetis_mark_parts(parts, msh, part)
211  else
212  call neko_error(rcode)
213  end if
214 
215  deallocate(part, xyz, vtxdist)
216 
217  end subroutine parmetis_partgeom
218 
220  subroutine parmetis_mark_parts(parts, msh, part)
221  type(mesh_fld_t), intent(inout) :: parts
222  type(mesh_t), intent(in) :: msh
223  integer(kind=M_INT), allocatable, intent(in) :: part(:)
224  integer :: i
225 
226  call mesh_field_init(parts, msh, 'partitions')
227 
228  do i = 1, msh%nelv
229  parts%data(i) = neko_idx(part(i))
230  end do
231 
232  end subroutine parmetis_mark_parts
233 
235  subroutine parmetis_wgt(msh, wgt, tpwgts, ubvec, nparts, ncon, weight)
236  type(mesh_t), intent(in) :: msh
237  integer(kind=M_INT), allocatable, intent(inout) :: wgt(:)
238  real(kind=m_real), allocatable, intent(inout) :: tpwgts(:)
239  real(kind=m_real), allocatable, intent(inout) :: ubvec(:)
240  integer, intent(in) :: nparts, ncon
241  type(mesh_fld_t), intent(in), optional :: weight
242  integer :: i
243 
244  if (present(weight)) then
245  do i = 1, msh%nelv
246  wgt(i) = parmetis_idx(weight%data(i))
247  end do
248  else
249  wgt = parmetis_idx(1)
250  end if
251 
252  do i = 1, (ncon * nparts)
253  tpwgts(i) = parmetis_real(1) / parmetis_real(nparts)
254  end do
255 
256  do i = 1, ncon
257  ubvec(i) = parmetis_real(1.05d0)
258  end do
259 
260  end subroutine parmetis_wgt
261 
263  subroutine parmetis_dist(dist, nelv)
264  integer(kind=M_INT), intent(inout) :: dist(0:pe_size)
265  integer, intent(in) :: nelv
266  integer(kind=M_INT) :: tmp, sum
267  integer :: i, ierr
268 
269  dist(pe_rank) = parmetis_idx(nelv)
270 
271  call mpi_allgather(nelv, 1, mpi_integer, dist, 1, &
272  mpi_integer, neko_comm, ierr)
273 
274  sum = dist(0)
275  dist(0) = 0
276  do i = 1, pe_size
277  tmp = dist(i)
278  dist(i) = sum
279  sum = tmp + sum
280  end do
281 
282  end subroutine parmetis_dist
283 
284 #else
285 
287  subroutine parmetis_partmeshkway(msh, parts, weights, nprts)
288  type(mesh_t), intent(inout) :: msh
289  type(mesh_fld_t), intent(inout) :: parts
290  type(mesh_fld_t), intent(in), optional :: weights
291  integer, intent(in), optional :: nprts
292  call neko_error('NEKO needs to be built with ParMETIS support')
293  end subroutine parmetis_partmeshkway
294 
297  subroutine parmetis_partgeom(msh, parts)
298  type(mesh_t), intent(inout) :: msh
299  type(mesh_fld_t), intent(inout) :: parts
300  call neko_error('NEKO needs to be built with ParMETIS support')
301  end subroutine parmetis_partgeom
302 
303 #endif
304 
305 end module parmetis
Definition: comm.F90:1
integer pe_rank
MPI rank.
Definition: comm.F90:28
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
integer pe_size
MPI size of communicator.
Definition: comm.F90:31
Defines a mesh field.
Definition: mesh_field.f90:35
subroutine, public mesh_field_init(fld, msh, fld_name)
Definition: mesh_field.f90:52
Defines a mesh.
Definition: mesh.f90:34
Interface to ParMETIS.
Definition: parmetis.F90:34
subroutine parmetis_mark_parts(parts, msh, part)
Fill mesh field according to new partitions.
Definition: parmetis.F90:221
integer, parameter, private metis_error
Definition: parmetis.F90:45
subroutine parmetis_dist(dist, nelv)
Compute the (parallel) vertex distribution of the dual graph.
Definition: parmetis.F90:264
subroutine, public parmetis_partmeshkway(msh, parts, weights, nprts)
Compute a k-way partitioning of a mesh msh.
Definition: parmetis.F90:111
integer, parameter, private metis_ok
Definition: parmetis.F90:45
subroutine parmetis_wgt(msh, wgt, tpwgts, ubvec, nparts, ncon, weight)
Setup weights and balance constraints for the dual graph.
Definition: parmetis.F90:236
subroutine, public parmetis_partgeom(msh, parts)
Compute a k-way partitioning of a mesh msh using a coordinated-based space-filing curves method.
Definition: parmetis.F90:182
Implements a point.
Definition: point.f90:35
Utilities.
Definition: utils.f90:35
A point in with coordinates .
Definition: point.f90:43