41 use,
intrinsic :: iso_c_binding
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
57 type(c_ptr),
value :: elmdist, eptr, eind, elmwgt, wgtflag, &
58 numflag, ncon, ncommonnodes, nparts, options, edgecut, part
60 type(c_ptr),
value :: tpwgts, ubvec
66 (vtxdist, ndims, xyz, part) &
67 bind(c, name=
'ParMETIS_V3_PartGeom_wrapper')
68 use,
intrinsic :: iso_c_binding
71 type(c_ptr),
value :: vtxdist, ndims, part
73 type(c_ptr),
value :: xyz
87 #ifdef HAVE_PARMETIS_REAL64
88 #define M_REAL c_double
89 #define parmetis_real(i) real((i), i8)
91 #define M_REAL c_float
92 #define parmetis_real(i) real((i), 4)
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)
99 #define M_INT c_int32_t
100 #define parmetis_idx(i) (i)
101 #define neko_idx(i) (i)
111 type(
mesh_t),
intent(inout) :: msh
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
127 ncommonnodes = 2**(msh%gdim - 1)
133 if (
present(nprts))
then
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))
145 if (
present(weights))
then
146 call parmetis_wgt(msh, elmwgt, tpwgts, ubvec, nparts, ncon, weights)
148 call parmetis_wgt(msh, elmwgt, tpwgts, ubvec, nparts, ncon)
153 eptr(i) = parmetis_idx(eptr(i - 1) + msh%npts)
159 eind(k) = parmetis_idx(msh%elements(i)%e%pts(j)%p%id() - 1)
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))
175 deallocate(elmdist, eptr, eind, part, elmwgt, tpwgts, ubvec)
182 type(
mesh_t),
intent(inout) :: msh
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
188 integer :: i, j, ierr, rcode
192 allocate(part(msh%nelv), xyz(ndims * 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))
207 c_loc(xyz), c_loc(part))
215 deallocate(part, xyz, vtxdist)
222 type(
mesh_t),
intent(in) :: msh
223 integer(kind=M_INT),
allocatable,
intent(in) :: part(:)
229 parts%data(i) = neko_idx(part(i))
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
244 if (
present(weight))
then
246 wgt(i) = parmetis_idx(weight%data(i))
249 wgt = parmetis_idx(1)
252 do i = 1, (ncon * nparts)
253 tpwgts(i) = parmetis_real(1) / parmetis_real(nparts)
257 ubvec(i) = parmetis_real(1.05d0)
264 integer(kind=M_INT),
intent(inout) :: dist(0:pe_size)
265 integer,
intent(in) :: nelv
266 integer(kind=M_INT) :: tmp, sum
269 dist(
pe_rank) = parmetis_idx(nelv)
271 call mpi_allgather(nelv, 1, mpi_integer, dist, 1, &
288 type(
mesh_t),
intent(inout) :: msh
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')
298 type(
mesh_t),
intent(inout) :: msh
300 call neko_error(
'NEKO needs to be built with ParMETIS support')
type(mpi_comm) neko_comm
MPI communicator.
integer pe_size
MPI size of communicator.
subroutine, public mesh_field_init(fld, msh, fld_name)
subroutine parmetis_mark_parts(parts, msh, part)
Fill mesh field according to new partitions.
integer, parameter, private metis_error
subroutine parmetis_dist(dist, nelv)
Compute the (parallel) vertex distribution of the dual graph.
subroutine, public parmetis_partmeshkway(msh, parts, weights, nprts)
Compute a k-way partitioning of a mesh msh.
integer, parameter, private metis_ok
subroutine parmetis_wgt(msh, wgt, tpwgts, ubvec, nparts, ncon, weight)
Setup weights and balance constraints for the dual graph.
subroutine, public parmetis_partgeom(msh, parts)
Compute a k-way partitioning of a mesh msh using a coordinated-based space-filing curves method.
A point in with coordinates .