13 #define ERROR(rval) if (0 .ne. rval) call exit(1)
20 #include "moab/MOABConfig.h"
23 # include "iMeshP_f.h"
33 integer numv, nume, nvpere
38 ibase_entityhandle,
pointer :: ents(:), verts(:)
39 ibase_entitysethandle root_set
40 TYPE(c_ptr) :: vertsptr, entsptr
42 real*8 coords(0:3*numv-1)
43 integer iconn(0:4*nume-1), gids(0:numv-1)
46 integer lgids(0:numv-1), lconn(0:4*nume-1)
47 real*8 lcoords(0:3*numv-1)
48 integer lnumv, lvids(0:numv-1), gvids(0:numv-1)
50 integer ic, ie, iv, istart, iend, ierr, indv, lnume, rank, sz
54 imeshp_partitionhandle imeshp
61 -1., -1., -1, 1., -1., -1., 1., 1., -1., -1., 1., -1., &
62 -1., -1., 1, 1., -1., 1., 1., 1., 1., -1., 1., 1. /
74 data ltp / imesh_quadrilateral /
84 call mpi_comm_size(mpi_comm_world, sz, ierr)
85 call mpi_comm_rank(mpi_comm_world, rank, ierr)
90 iend = istart + lnume - 1
91 if (rank .eq. sz-1)
then
93 lnume = iend - istart + 1
104 indv = iconn(lvpe*ie + iv)
105 if (lvids(indv) .eq. -1)
then
108 lcoords(3*lnumv+ic) = coords(3*indv+ic)
111 gvids(lnumv) = 1+indv
113 lconn(lvpe*(ie-istart)+iv) = lvids(indv)
122 call create_mesh(imesh, imeshp, mpi_comm_world, lnumv, lnume, gvids, lvpe, ltp, lcoords, lconn, &
123 vertsptr, entsptr, ierr)
125 call c_f_pointer(vertsptr, verts, [lnumv])
126 call c_f_pointer(entsptr, ents, [lnume])
129 call imesh_getrootset(%VAL(imesh), root_set, ierr)
133 call imeshp_getnumoftypeall(%VAL(imesh), %VAL(imeshp), %VAL(root_set), %VAL(ibase_vertex), iv, ierr)
135 call imeshp_getnumoftypeall(%VAL(imesh), %VAL(imeshp), %VAL(root_set), %VAL(ibase_face), ie, ierr)
137 if (rank .eq. 0)
then
138 write(0,*)
"Number of vertices = ", iv
139 write(0,*)
"Number of entities = ", ie
144 call mpi_finalize(ierr)
146 write(0, *)
"compile with MPI for better experience\n"
152 subroutine create_mesh( &
156 comm, numv, nume, vgids, nvpe, tp, posn, iconn, &
158 vertsPtr, entsPtr, ierr)
171 # include "iMeshP_f.h"
176 TYPE(c_ptr) :: vertsptr, entsptr
177 integer numv, nume, nvpe, vgids(0:*), iconn(0:*), ierr, tp
180 imeshp_partitionhandle imeshp
185 integer comm_sz, comm_rank, numa, numo, iv, ie
186 TYPE(c_ptr) :: statsptr
187 integer,
allocatable,
target :: stats(:)
190 ibase_entityhandle,
pointer :: verts(:), ents(:)
191 ibase_entityhandle,
allocatable :: conn(:)
192 ibase_entitysethandle root_set
193 ibase_entitysethandle file_set
195 ibase_handle_t mpi_comm_c
196 TYPE(c_ptr) :: partsptr
197 imeshp_parthandle,
pointer :: parts(:)
198 imeshp_parthandle part
199 integer partsa, partso
200 character (len=10) filename
204 if (imesh .eq. 0)
then
205 call imesh_newmesh(
"MOAB", imesh, ierr)
209 if (imeshp .eq. 0)
then
210 call imeshp_getcommunicator(%VAL(imesh), mpi_comm_world, mpi_comm_c, ierr)
212 call imeshp_createpartitionall(%VAL(imesh), %VAL(mpi_comm_c), imeshp, ierr)
214 call imeshp_createpart(%VAL(imesh), %VAL(imeshp), part, ierr)
218 call imeshp_getlocalparts(%VAL(imesh), %VAL(imeshp), partsptr, partsa, partso, ierr)
220 call c_f_pointer(partsptr, parts, [partso])
223 call mpi_comm_rank(comm, comm_rank, ierr)
225 call mpi_comm_size(comm, comm_sz, ierr)
231 call imesh_createvtxarr(%VAL(imesh), %VAL(numv), %VAL(ibase_interleaved), posn, %VAL(3*numv), &
232 vertsptr, numa, numo, ierr)
236 allocate (conn(0:nvpe*nume-1))
237 call c_f_pointer(vertsptr, verts, [numv])
238 do i = 0, nvpe*nume-1
239 conn(i) = verts(1+iconn(i))
243 allocate(stats(0:nume-1))
244 statsptr = c_loc(stats(0))
245 call imesh_createentarr(%VAL(imesh), %VAL(tp), conn, %VAL(nvpe*nume), &
246 entsptr, numa, numo, statsptr, numa, numo, ierr)
254 call c_f_pointer(entsptr, ents, [numo])
255 call imesh_addentarrtoset(%VAL(imesh), ents, %VAL(numo), %VAL(part), ierr)
258 call imesh_gettaghandle(%VAL(imesh),
"GLOBAL_ID", tagh, ierr, %VAL(9))
259 if (ibase_success .ne. ierr)
then
261 call imesh_createtag(%VAL(imesh),
"GLOBAL_ID", %VAL(ibase_integer), tagh, ierr)
264 call imesh_setintarrdata(%VAL(imesh), verts, %VAL(numv), %VAL(tagh), vgids, %VAL(numv), ierr)
267 call imeshp_syncmeshall(%VAL(imesh), %VAL(imeshp), ierr)
270 call imeshp_createghostentsall(%VAL(imesh), %VAL(imeshp), %VAL(2), %VAL(1), %VAL(1), %VAL(0), ierr)
276 end subroutine create_mesh