xref: /petsc/src/binding/petsc4py/test/test_dmplex.py (revision 6220ef827cfb29f143b71dc8d0e16ee8b5f4575d)
1import petsc4py
2from petsc4py import PETSc
3import unittest
4import os
5import filecmp
6import numpy as np
7import importlib
8
9# --------------------------------------------------------------------
10
11ERR_ARG_OUTOFRANGE = 63
12
13
14class BaseTestPlex:
15    COMM = PETSc.COMM_WORLD
16    DIM = 1
17    CELLS = [[0, 1], [1, 2]]
18    COORDS = [[0.0], [0.5], [1.0]]
19    COMP = 1
20    DOFS = [1, 0]
21
22    def setUp(self):
23        self.plex = PETSc.DMPlex().createFromCellList(
24            self.DIM, self.CELLS, self.COORDS, comm=self.COMM
25        )
26
27    def tearDown(self):
28        self.plex.destroy()
29        self.plex = None
30        PETSc.garbage_cleanup()
31
32    def testTopology(self):
33        rank = self.COMM.rank
34        dim = self.plex.getDimension()
35        pStart, pEnd = self.plex.getChart()
36        cStart, cEnd = self.plex.getHeightStratum(0)
37        vStart, vEnd = self.plex.getDepthStratum(0)
38        numDepths = self.plex.getLabelSize('depth')
39        coords_raw = self.plex.getCoordinates().getArray()
40        coords = np.reshape(coords_raw, (vEnd - vStart, dim))
41        self.assertEqual(dim, self.DIM)
42        self.assertEqual(numDepths, self.DIM + 1)
43        if rank == 0 and self.CELLS is not None:
44            self.assertEqual(cEnd - cStart, len(self.CELLS))
45        if rank == 0 and self.COORDS is not None:
46            self.assertEqual(vEnd - vStart, len(self.COORDS))
47            self.assertTrue((coords == self.COORDS).all())
48
49    def testClosure(self):
50        pStart, pEnd = self.plex.getChart()
51        for p in range(pStart, pEnd):
52            closure = self.plex.getTransitiveClosure(p)[0]
53            for c in closure:
54                cone = self.plex.getCone(c)
55                self.assertEqual(self.plex.getConeSize(c), len(cone))
56                for i in cone:
57                    self.assertIn(i, closure)
58            star = self.plex.getTransitiveClosure(p, useCone=False)[0]
59            for s in star:
60                support = self.plex.getSupport(s)
61                self.assertEqual(self.plex.getSupportSize(s), len(support))
62                for i in support:
63                    self.assertIn(i, star)
64
65    def testAdjacency(self):
66        PETSc.DMPlex.setAdjacencyUseAnchors(self.plex, False)
67        flag = PETSc.DMPlex.getAdjacencyUseAnchors(self.plex)
68        self.assertFalse(flag)
69        PETSc.DMPlex.setAdjacencyUseAnchors(self.plex, True)
70        flag = PETSc.DMPlex.getAdjacencyUseAnchors(self.plex)
71        self.assertTrue(flag)
72        PETSc.DMPlex.setBasicAdjacency(self.plex, False, False)
73        flagA, flagB = PETSc.DMPlex.getBasicAdjacency(self.plex)
74        self.assertFalse(flagA)
75        self.assertFalse(flagB)
76        PETSc.DMPlex.setBasicAdjacency(self.plex, True, True)
77        flagA, flagB = PETSc.DMPlex.getBasicAdjacency(self.plex)
78        self.assertTrue(flagA)
79        self.assertTrue(flagB)
80        pStart, pEnd = self.plex.getChart()
81        for p in range(pStart, pEnd):
82            adjacency = self.plex.getAdjacency(p)
83            self.assertTrue(p in adjacency)
84            self.assertTrue(len(adjacency) > 1)
85
86    def testSectionDofs(self):
87        self.plex.setNumFields(1)
88        section = self.plex.createSection([self.COMP], [self.DOFS])
89        size = section.getStorageSize()
90        entity_dofs = [
91            self.plex.getStratumSize('depth', d) * self.DOFS[d]
92            for d in range(self.DIM + 1)
93        ]
94        self.assertEqual(sum(entity_dofs), size)
95
96    def testSectionClosure(self):
97        section = self.plex.createSection([self.COMP], [self.DOFS])
98        self.plex.setSection(section)
99        vec = self.plex.createLocalVec()
100        pStart, pEnd = self.plex.getChart()
101        for p in range(pStart, pEnd):
102            for i in range(section.getDof(p)):
103                off = section.getOffset(p)
104                vec.setValue(off + i, p)
105
106        for p in range(pStart, pEnd):
107            point_closure = self.plex.getTransitiveClosure(p)[0]
108            dof_closure = self.plex.vecGetClosure(section, vec, p)
109            for p in dof_closure:
110                self.assertIn(p, point_closure)
111
112    def testBoundaryLabel(self):
113        pStart, pEnd = self.plex.getChart()
114        if pEnd - pStart == 0:
115            return
116
117        self.assertFalse(self.plex.hasLabel('boundary'))
118        self.plex.markBoundaryFaces('boundary')
119        self.assertTrue(self.plex.hasLabel('boundary'))
120
121        faces = self.plex.getStratumIS('boundary', 1)
122        for f in faces.getIndices():
123            points, orient = self.plex.getTransitiveClosure(f, useCone=True)
124            for p in points:
125                self.plex.setLabelValue('boundary', p, 1)
126
127        for p in range(pStart, pEnd):
128            if self.plex.getLabelValue('boundary', p) != 1:
129                self.plex.setLabelValue('boundary', p, 2)
130
131        numBoundary = self.plex.getStratumSize('boundary', 1)
132        numInterior = self.plex.getStratumSize('boundary', 2)
133        self.assertNotEqual(numBoundary, pEnd - pStart)
134        self.assertNotEqual(numInterior, pEnd - pStart)
135        self.assertEqual(numBoundary + numInterior, pEnd - pStart)
136
137    def testMetric(self):
138        if self.DIM == 1:
139            return
140        self.plex.distribute()
141        if self.CELLS is None and not self.plex.isSimplex():
142            return
143        self.plex.orient()
144
145        h_min = 1.0e-30
146        h_max = 1.0e30
147        a_max = 1.0e10
148        target = 8.0
149        p = 1.0
150        beta = 1.3
151        hausd = 0.01
152        self.plex.metricSetUniform(False)
153        self.plex.metricSetIsotropic(False)
154        self.plex.metricSetRestrictAnisotropyFirst(False)
155        self.plex.metricSetNoInsertion(False)
156        self.plex.metricSetNoSwapping(False)
157        self.plex.metricSetNoMovement(False)
158        self.plex.metricSetNoSurf(False)
159        self.plex.metricSetVerbosity(-1)
160        self.plex.metricSetNumIterations(3)
161        self.plex.metricSetMinimumMagnitude(h_min)
162        self.plex.metricSetMaximumMagnitude(h_max)
163        self.plex.metricSetMaximumAnisotropy(a_max)
164        self.plex.metricSetTargetComplexity(target)
165        self.plex.metricSetNormalizationOrder(p)
166        self.plex.metricSetGradationFactor(beta)
167        self.plex.metricSetHausdorffNumber(hausd)
168
169        self.assertFalse(self.plex.metricIsUniform())
170        self.assertFalse(self.plex.metricIsIsotropic())
171        self.assertFalse(self.plex.metricRestrictAnisotropyFirst())
172        self.assertFalse(self.plex.metricNoInsertion())
173        self.assertFalse(self.plex.metricNoSwapping())
174        self.assertFalse(self.plex.metricNoMovement())
175        self.assertFalse(self.plex.metricNoSurf())
176        self.assertTrue(self.plex.metricGetVerbosity() == -1)
177        self.assertTrue(self.plex.metricGetNumIterations() == 3)
178        self.assertTrue(np.isclose(self.plex.metricGetMinimumMagnitude(), h_min))
179        self.assertTrue(np.isclose(self.plex.metricGetMaximumMagnitude(), h_max))
180        self.assertTrue(np.isclose(self.plex.metricGetMaximumAnisotropy(), a_max))
181        self.assertTrue(np.isclose(self.plex.metricGetTargetComplexity(), target))
182        self.assertTrue(np.isclose(self.plex.metricGetNormalizationOrder(), p))
183        self.assertTrue(np.isclose(self.plex.metricGetGradationFactor(), beta))
184        self.assertTrue(np.isclose(self.plex.metricGetHausdorffNumber(), hausd))
185
186        metric1 = self.plex.metricCreateUniform(0.5)
187        metric2 = self.plex.metricCreateUniform(1.0)
188        metric = self.plex.metricCreate()
189        det = self.plex.metricDeterminantCreate()
190        self.plex.metricAverage2(metric1, metric2, metric)
191        metric1.array[:] *= 1.5
192        self.assertTrue(np.allclose(metric.array, metric1.array))
193        self.plex.metricIntersection2(metric1, metric2, metric)
194        self.assertTrue(np.allclose(metric.array, metric2.array))
195        self.plex.metricEnforceSPD(metric, metric1, det[0])
196        self.assertTrue(np.allclose(metric.array, metric1.array))
197
198        if self.DIM == 2 and PETSc.COMM_WORLD.getSize() > 6:
199            # Error with 7 processes in 2D: normalization factor is -1
200            return
201
202        self.plex.metricNormalize(
203            metric, metric1, det[0], restrictSizes=False, restrictAnisotropy=False
204        )
205        metric2.scale(pow(target, 2.0 / self.DIM))
206        self.assertTrue(np.allclose(metric1.array, metric2.array))
207
208    def testAdapt(self):
209        if self.DIM == 1:
210            return
211        if self.DIM == 3 and PETSc.COMM_WORLD.getSize() > 4:
212            # Error with 5 processes in 3D
213            # ----------------------------
214            # Warning: MMG5_mmgIntextmet: Unable to diagonalize at least 1 metric.
215            # Error: MMG3D_defsiz_ani: unable to intersect metrics at point 8.
216            # Metric undefined. Exit program.
217            # MMG remeshing problem. Exit program.
218            return
219        self.plex.orient()
220        plex = self.plex.refine()
221        plex.distribute()
222        if self.CELLS is None and not plex.isSimplex():
223            return
224        if sum(self.DOFS) > 1:
225            return
226        metric = plex.metricCreateUniform(9.0)
227        try:
228            newplex = plex.adaptMetric(metric, '')
229            plex.destroy()
230            newplex.destroy()
231        except PETSc.Error as exc:
232            plex.destroy()
233            if exc.ierr != ERR_ARG_OUTOFRANGE:
234                raise
235
236    def testNatural(self):
237        dim = self.plex.getDimension()
238        ct = self.plex.getCellType(0)
239        fe = PETSc.FE().createByCell(dim, 1, ct)
240        self.plex.setField(0, fe)
241        self.plex.createDS()
242        self.plex.setUseNatural(True)
243        self.plex.distribute()
244        self.plex.view()
245        gv = self.plex.createGlobalVec()
246        nv = self.plex.createNaturalVec()
247        self.plex.globalToNaturalBegin(gv, nv)
248        self.plex.globalToNaturalEnd(gv, nv)
249        self.plex.naturalToGlobalBegin(nv, gv)
250        self.plex.naturalToGlobalEnd(nv, gv)
251
252# --------------------------------------------------------------------
253
254
255class BaseTestPlex_2D(BaseTestPlex):
256    DIM = 2
257    CELLS = [
258        [0, 1, 3],
259        [1, 3, 4],
260        [1, 2, 4],
261        [2, 4, 5],
262        [3, 4, 6],
263        [4, 6, 7],
264        [4, 5, 7],
265        [5, 7, 8],
266    ]
267    COORDS = [
268        [0.0, 0.0],
269        [0.5, 0.0],
270        [1.0, 0.0],
271        [0.0, 0.5],
272        [0.5, 0.5],
273        [1.0, 0.5],
274        [0.0, 1.0],
275        [0.5, 1.0],
276        [1.0, 1.0],
277    ]
278    DOFS = [1, 0, 0]
279
280
281class BaseTestPlex_3D(BaseTestPlex):
282    DIM = 3
283    CELLS = [
284        [0, 2, 3, 7],
285        [0, 2, 6, 7],
286        [0, 4, 6, 7],
287        [0, 1, 3, 7],
288        [0, 1, 5, 7],
289        [0, 4, 5, 7],
290    ]
291    COORDS = [
292        [0.0, 0.0, 0.0],
293        [1.0, 0.0, 0.0],
294        [0.0, 1.0, 0.0],
295        [1.0, 1.0, 0.0],
296        [0.0, 0.0, 1.0],
297        [1.0, 0.0, 1.0],
298        [0.0, 1.0, 1.0],
299        [1.0, 1.0, 1.0],
300    ]
301    DOFS = [1, 0, 0, 0]
302
303
304# --------------------------------------------------------------------
305
306
307class TestPlex_1D(BaseTestPlex, unittest.TestCase):
308    pass
309
310
311class TestPlex_2D(BaseTestPlex_2D, unittest.TestCase):
312    def testTransform(self):
313        plex = self.plex
314        cstart, cend = plex.getHeightStratum(0)
315        tr = PETSc.DMPlexTransform().create(comm=PETSc.COMM_WORLD)
316        tr.setType(PETSc.DMPlexTransformType.REFINEALFELD)
317        tr.setDM(plex)
318        tr.setUp()
319        newplex = tr.apply(plex)
320        tr.destroy()
321        newcstart, newcend = newplex.getHeightStratum(0)
322        newplex.destroy()
323        self.assertTrue((newcend - newcstart) == 3 * (cend - cstart))
324
325
326class TestPlex_3D(BaseTestPlex_3D, unittest.TestCase):
327    pass
328
329
330class TestPlex_2D_P3(BaseTestPlex_2D, unittest.TestCase):
331    DOFS = [1, 2, 1]
332
333
334class TestPlex_3D_P3(BaseTestPlex_3D, unittest.TestCase):
335    DOFS = [1, 2, 1, 0]
336
337
338class TestPlex_3D_P4(BaseTestPlex_3D, unittest.TestCase):
339    DOFS = [1, 3, 3, 1]
340
341
342class TestPlex_2D_BoxTensor(BaseTestPlex_2D, unittest.TestCase):
343    CELLS = None
344    COORDS = None
345
346    def setUp(self):
347        self.plex = PETSc.DMPlex().createBoxMesh([3, 3], simplex=False)
348
349
350class TestPlex_3D_BoxTensor(BaseTestPlex_3D, unittest.TestCase):
351    CELLS = None
352    COORDS = None
353
354    def setUp(self):
355        self.plex = PETSc.DMPlex().createBoxMesh([3, 3, 3], simplex=False)
356
357
358# FIXME
359try:
360    raise PETSc.Error
361    PETSc.DMPlex().createBoxMesh([2, 2], simplex=True, comm=PETSc.COMM_SELF).destroy()
362except PETSc.Error:
363    pass
364else:
365
366    class TestPlex_2D_Box(BaseTestPlex_2D, unittest.TestCase):
367        CELLS = None
368        COORDS = None
369
370        def setUp(self):
371            self.plex = PETSc.DMPlex().createBoxMesh([1, 1], simplex=True)
372
373    class TestPlex_2D_Boundary(BaseTestPlex_2D, unittest.TestCase):
374        CELLS = None
375        COORDS = None
376
377        def setUp(self):
378            boundary = PETSc.DMPlex().create(self.COMM)
379            boundary.createSquareBoundary([0.0, 0.0], [1.0, 1.0], [2, 2])
380            boundary.setDimension(self.DIM - 1)
381            self.plex = PETSc.DMPlex().generate(boundary)
382
383    class TestPlex_3D_Box(BaseTestPlex_3D, unittest.TestCase):
384        CELLS = None
385        COORDS = None
386
387        def setUp(self):
388            self.plex = PETSc.DMPlex().createBoxMesh([1, 1, 1], simplex=True)
389
390    class TestPlex_3D_Boundary(BaseTestPlex_3D, unittest.TestCase):
391        CELLS = None
392        COORDS = None
393
394        def setUp(self):
395            boundary = PETSc.DMPlex().create(self.COMM)
396            boundary.createCubeBoundary([0.0, 0.0, 0.0], [1.0, 1.0, 1.0], [1, 1, 1])
397            boundary.setDimension(self.DIM - 1)
398            self.plex = PETSc.DMPlex().generate(boundary)
399
400# --------------------------------------------------------------------
401
402PETSC_DIR = petsc4py.get_config()['PETSC_DIR']
403
404
405def check_dtype(method):
406    def wrapper(self, *args, **kwargs):
407        if PETSc.ScalarType is PETSc.ComplexType:
408            return None
409        return method(self, *args, **kwargs)
410
411    return wrapper
412
413
414def check_package(method):
415    def wrapper(self, *args, **kwargs):
416        if not PETSc.Sys.hasExternalPackage('hdf5'):
417            return None
418        if self.PARTITIONERTYPE != 'simple' and not PETSc.Sys.hasExternalPackage(
419            self.PARTITIONERTYPE
420        ):
421            return None
422        return method(self, *args, **kwargs)
423
424    return wrapper
425
426
427def check_nsize(method):
428    def wrapper(self, *args, **kwargs):
429        if PETSc.COMM_WORLD.size != self.NSIZE:
430            return None
431        return method(self, *args, **kwargs)
432
433    return wrapper
434
435
436class BaseTestPlexHDF5:
437    NSIZE = 4
438    NTIMES = 3
439
440    def tearDown(self):
441        if not PETSc.COMM_WORLD.rank:
442            if os.path.exists(self.outfile()):
443                os.remove(self.outfile())
444            if os.path.exists(self.tmp_output_file()):
445                os.remove(self.tmp_output_file())
446
447    def _name(self):
448        return f'{self.SUFFIX}_outformat-{self.OUTFORMAT}_{self.PARTITIONERTYPE}'
449
450    def infile(self):
451        return os.path.join(
452            PETSC_DIR, 'share/petsc/datafiles/', 'meshes/blockcylinder-50.h5'
453        )
454
455    def outfile(self):
456        return os.path.join('./temp_test_dmplex_%s.h5' % self._name())
457
458    def informat(self):
459        return PETSc.Viewer.Format.HDF5_XDMF
460
461    def outformat(self):
462        d = {
463            'hdf5_petsc': PETSc.Viewer.Format.HDF5_PETSC,
464            'hdf5_xdmf': PETSc.Viewer.Format.HDF5_XDMF,
465        }
466        return d[self.OUTFORMAT]
467
468    def partitionerType(self):
469        d = {
470            'simple': PETSc.Partitioner.Type.SIMPLE,
471            'ptscotch': PETSc.Partitioner.Type.PTSCOTCH,
472            'parmetis': PETSc.Partitioner.Type.PARMETIS,
473        }
474        return d[self.PARTITIONERTYPE]
475
476    def ref_output_file(self):
477        return os.path.join(
478            PETSC_DIR,
479            'src/dm/impls/plex/tutorials/',
480            'output/ex5_%s.out' % self._name(),
481        )
482
483    def tmp_output_file(self):
484        return os.path.join('./temp_test_dmplex_%s.out' % self._name())
485
486    def outputText(self, msg, comm):
487        if not comm.rank:
488            with open(self.tmp_output_file(), 'a') as f:
489                f.write(msg)
490
491    def outputPlex(self, plex):
492        txtvwr = PETSc.Viewer().createASCII(
493            self.tmp_output_file(), mode='a', comm=plex.comm
494        )
495        plex.view(viewer=txtvwr)
496        txtvwr.destroy()
497
498    @check_dtype
499    @check_package
500    @check_nsize
501    def testViewLoadCycle(self):
502        if importlib.util.find_spec('mpi4py') is None:
503            self.skipTest('mpi4py')  # throws special exception to signal test skip
504        grank = PETSc.COMM_WORLD.rank
505        for i in range(self.NTIMES):
506            if i == 0:
507                infname = self.infile()
508                informt = self.informat()
509            else:
510                infname = self.outfile()
511                informt = self.outformat()
512            if self.HETEROGENEOUS:
513                mycolor = grank > self.NTIMES - i
514            else:
515                mycolor = 0
516            mpicomm = PETSc.COMM_WORLD.tompi4py()
517            comm = PETSc.Comm(comm=mpicomm.Split(color=mycolor, key=grank))
518            if mycolor == 0:
519                self.outputText('Begin cycle %d\n' % i, comm)
520                plex = PETSc.DMPlex()
521                vwr = PETSc.ViewerHDF5()
522                # Create plex
523                plex.create(comm=comm)
524                plex.setName('DMPlex Object')
525                # Load data from XDMF into dm in parallel
526                vwr.create(infname, mode='r', comm=comm)
527                vwr.pushFormat(format=informt)
528                plex.load(viewer=vwr)
529                plex.setOptionsPrefix('loaded_')
530                plex.distributeSetDefault(False)
531                plex.setFromOptions()
532                vwr.popFormat()
533                vwr.destroy()
534                self.outputPlex(plex)
535                # Test DM is indeed distributed
536                flg = plex.isDistributed()
537                self.outputText(
538                    'Loaded mesh distributed? %s\n' % str(flg).upper(), comm
539                )
540                # Interpolate
541                plex.interpolate()
542                plex.setOptionsPrefix('interpolated_')
543                plex.setFromOptions()
544                self.outputPlex(plex)
545                # Redistribute
546                part = plex.getPartitioner()
547                part.setType(self.partitionerType())
548                sf = plex.distribute(overlap=0)
549                if sf:
550                    sf.destroy()
551                part.destroy()
552                plex.setName('DMPlex Object')
553                plex.setOptionsPrefix('redistributed_')
554                plex.setFromOptions()
555                self.outputPlex(plex)
556                # Save redistributed dm to XDMF in parallel
557                vwr.create(self.outfile(), mode='w', comm=comm)
558                vwr.pushFormat(format=self.outformat())
559                plex.setName('DMPlex Object')
560                plex.view(viewer=vwr)
561                vwr.popFormat()
562                vwr.destroy()
563                # Destroy plex
564                plex.destroy()
565                self.outputText('End   cycle %d\n--------\n' % i, comm)
566            comm.tompi4py().Free()
567            PETSc.COMM_WORLD.Barrier()
568        # Check that the output is identical to that of plex/tutorial/ex5.c.
569        self.assertTrue(
570            filecmp.cmp(self.tmp_output_file(), self.ref_output_file(), shallow=False),
571            f'Contents of the files not the same. Reference file: {self.ref_output_file()}',
572        )
573        PETSc.COMM_WORLD.Barrier()
574
575
576class BaseTestPlexHDF5Homogeneous(BaseTestPlexHDF5):
577    """Test save on N / load on N."""
578
579    SUFFIX = 0
580    HETEROGENEOUS = False
581
582
583class BaseTestPlexHDF5Heterogeneous(BaseTestPlexHDF5):
584    """Test save on N / load on M."""
585
586    SUFFIX = 1
587    HETEROGENEOUS = True
588
589
590class TestPlexHDF5PETSCSimpleHomogeneous(
591    BaseTestPlexHDF5Homogeneous, unittest.TestCase
592):
593    OUTFORMAT = 'hdf5_petsc'
594    PARTITIONERTYPE = 'simple'
595
596
597"""
598Skipping. PTScotch produces different distributions when run
599in a sequence in a single session.
600
601class TestPlexHDF5PETSCPTScotchHomogeneous(BaseTestPlexHDF5Homogeneous,
602                                           unittest.TestCase):
603    OUTFORMAT = "hdf5_petsc"
604    PARTITIONERTYPE = "ptscotch"
605"""
606
607
608class TestPlexHDF5PETSCParmetisHomogeneous(
609    BaseTestPlexHDF5Homogeneous, unittest.TestCase
610):
611    OUTFORMAT = 'hdf5_petsc'
612    PARTITIONERTYPE = 'parmetis'
613
614
615class TestPlexHDF5XDMFSimpleHomogeneous(BaseTestPlexHDF5Homogeneous, unittest.TestCase):
616    OUTFORMAT = 'hdf5_xdmf'
617    PARTITIONERTYPE = 'simple'
618
619
620"""
621Skipping. PTScotch produces different distributions when run
622in a sequence in a single session.
623
624class TestPlexHDF5XDMFPTScotchHomogeneous(BaseTestPlexHDF5Homogeneous,
625                                          unittest.TestCase):
626    OUTFORMAT = "hdf5_xdmf"
627    PARTITIONERTYPE = "ptscotch"
628"""
629
630
631class TestPlexHDF5XDMFParmetisHomogeneous(
632    BaseTestPlexHDF5Homogeneous, unittest.TestCase
633):
634    OUTFORMAT = 'hdf5_xdmf'
635    PARTITIONERTYPE = 'parmetis'
636
637
638class TestPlexHDF5PETSCSimpleHeterogeneous(
639    BaseTestPlexHDF5Heterogeneous, unittest.TestCase
640):
641    OUTFORMAT = 'hdf5_petsc'
642    PARTITIONERTYPE = 'simple'
643
644
645"""
646Skipping. PTScotch produces different distributions when run
647in a sequence in a single session.
648
649class TestPlexHDF5PETSCPTScotchHeterogeneous(BaseTestPlexHDF5Heterogeneous,
650                                             unittest.TestCase):
651    OUTFORMAT = "hdf5_petsc"
652    PARTITIONERTYPE = "ptscotch"
653"""
654
655
656class TestPlexHDF5PETSCParmetisHeterogeneous(
657    BaseTestPlexHDF5Heterogeneous, unittest.TestCase
658):
659    OUTFORMAT = 'hdf5_petsc'
660    PARTITIONERTYPE = 'parmetis'
661
662
663class TestPlexHDF5XDMFSimpleHeterogeneous(
664    BaseTestPlexHDF5Heterogeneous, unittest.TestCase
665):
666    OUTFORMAT = 'hdf5_xdmf'
667    PARTITIONERTYPE = 'simple'
668
669
670class TestPlexHDF5XDMFPTScotchHeterogeneous(
671    BaseTestPlexHDF5Heterogeneous, unittest.TestCase
672):
673    OUTFORMAT = 'hdf5_xdmf'
674    PARTITIONERTYPE = 'ptscotch'
675
676
677class TestPlexHDF5XDMFParmetisHeterogeneous(
678    BaseTestPlexHDF5Heterogeneous, unittest.TestCase
679):
680    OUTFORMAT = 'hdf5_xdmf'
681    PARTITIONERTYPE = 'parmetis'
682
683
684# --------------------------------------------------------------------
685
686if __name__ == '__main__':
687    unittest.main()
688