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