xref: /petsc/src/binding/petsc4py/test/test_dmplex.py (revision d6a8d7c6470df1f76c26d5e390a3bca0663924ae)
1import petsc4py
2from petsc4py import PETSc
3import unittest
4import os
5import filecmp
6import numpy as np
7
8# --------------------------------------------------------------------
9
10ERR_SUP = 56
11
12class BaseTestPlex(object):
13
14    COMM = PETSc.COMM_WORLD
15    DIM = 1
16    CELLS = [[0, 1], [1, 2]]
17    COORDS = [[0.], [0.5], [1.]]
18    COMP = 1
19    DOFS = [1, 0]
20
21    def setUp(self):
22        self.plex = PETSc.DMPlex().createFromCellList(self.DIM,
23                                                      self.CELLS,
24                                                      self.COORDS,
25                                                      comm=self.COMM)
26
27    def tearDown(self):
28        self.plex.destroy()
29        self.plex = None
30
31    def testTopology(self):
32        rank = self.COMM.rank
33        dim = self.plex.getDimension()
34        pStart, pEnd = self.plex.getChart()
35        cStart, cEnd = self.plex.getHeightStratum(0)
36        vStart, vEnd = self.plex.getDepthStratum(0)
37        numDepths = self.plex.getLabelSize("depth")
38        coords_raw = self.plex.getCoordinates().getArray()
39        coords = np.reshape(coords_raw, (vEnd - vStart, dim))
40        self.assertEqual(dim, self.DIM)
41        self.assertEqual(numDepths, self.DIM+1)
42        if rank == 0 and self.CELLS is not None:
43            self.assertEqual(cEnd-cStart, len(self.CELLS))
44        if rank == 0 and self.COORDS is not None:
45            self.assertEqual(vEnd-vStart, len(self.COORDS))
46            self.assertTrue((coords == self.COORDS).all())
47
48    def testClosure(self):
49        pStart, pEnd = self.plex.getChart()
50        for p in range(pStart, pEnd):
51            closure = self.plex.getTransitiveClosure(p)[0]
52            for c in closure:
53                cone = self.plex.getCone(c)
54                self.assertEqual(self.plex.getConeSize(c), len(cone))
55                for i in cone:
56                    self.assertIn(i, closure)
57            star = self.plex.getTransitiveClosure(p, useCone=False)[0]
58            for s in star:
59                support = self.plex.getSupport(s)
60                self.assertEqual(self.plex.getSupportSize(s), len(support))
61                for i in support:
62                    self.assertIn(i, star)
63
64    def testAdjacency(self):
65        PETSc.DMPlex.setAdjacencyUseAnchors(self.plex, False)
66        flag = PETSc.DMPlex.getAdjacencyUseAnchors(self.plex)
67        self.assertFalse(flag)
68        PETSc.DMPlex.setAdjacencyUseAnchors(self.plex, True)
69        flag = PETSc.DMPlex.getAdjacencyUseAnchors(self.plex)
70        self.assertTrue(flag)
71        PETSc.DMPlex.setBasicAdjacency(self.plex, False, False)
72        flagA, flagB = PETSc.DMPlex.getBasicAdjacency(self.plex)
73        self.assertFalse(flagA)
74        self.assertFalse(flagB)
75        PETSc.DMPlex.setBasicAdjacency(self.plex, True, True)
76        flagA, flagB = PETSc.DMPlex.getBasicAdjacency(self.plex)
77        self.assertTrue(flagA)
78        self.assertTrue(flagB)
79        pStart, pEnd = self.plex.getChart()
80        for p in range(pStart, pEnd):
81            adjacency = self.plex.getAdjacency(p)
82            self.assertTrue(p in adjacency)
83            self.assertTrue(len(adjacency) > 1)
84
85    def testSectionDofs(self):
86        self.plex.setNumFields(1)
87        section = self.plex.createSection([self.COMP], [self.DOFS])
88        size = section.getStorageSize()
89        entity_dofs = [self.plex.getStratumSize("depth", d) *
90                       self.DOFS[d] for d in range(self.DIM+1)]
91        self.assertEqual(sum(entity_dofs), size)
92
93    def testSectionClosure(self):
94        section = self.plex.createSection([self.COMP], [self.DOFS])
95        self.plex.setSection(section)
96        vec = self.plex.createLocalVec()
97        pStart, pEnd = self.plex.getChart()
98        for p in range(pStart, pEnd):
99            for i in range(section.getDof(p)):
100                off = section.getOffset(p)
101                vec.setValue(off+i, p)
102
103        for p in range(pStart, pEnd):
104            point_closure = self.plex.getTransitiveClosure(p)[0]
105            dof_closure = self.plex.vecGetClosure(section, vec, p)
106            for p in dof_closure:
107                self.assertIn(p, point_closure)
108
109    def testBoundaryLabel(self):
110        pStart, pEnd = self.plex.getChart()
111        if (pEnd - pStart == 0): return
112
113        self.assertFalse(self.plex.hasLabel("boundary"))
114        self.plex.markBoundaryFaces("boundary")
115        self.assertTrue(self.plex.hasLabel("boundary"))
116
117        faces = self.plex.getStratumIS("boundary", 1)
118        for f in faces.getIndices():
119            points, orient = self.plex.getTransitiveClosure(f, useCone=True)
120            for p in points:
121                self.plex.setLabelValue("boundary", p, 1)
122
123        for p in range(pStart, pEnd):
124            if self.plex.getLabelValue("boundary", p) != 1:
125                self.plex.setLabelValue("boundary", p, 2)
126
127        numBoundary = self.plex.getStratumSize("boundary", 1)
128        numInterior = self.plex.getStratumSize("boundary", 2)
129        self.assertNotEqual(numBoundary, pEnd - pStart)
130        self.assertNotEqual(numInterior, pEnd - pStart)
131        self.assertEqual(numBoundary + numInterior, pEnd - pStart)
132
133    def testMetric(self):
134        if self.DIM == 1: return
135        self.plex.distribute()
136        if self.CELLS is None and not self.plex.isSimplex(): return
137
138        h_min = 1.0e-30
139        h_max = 1.0e+30
140        a_max = 1.0e+10
141        target = 10.0
142        p = 1.0
143        self.plex.metricSetIsotropic(False)
144        self.plex.metricSetRestrictAnisotropyFirst(False)
145        self.plex.metricSetMinimumMagnitude(h_min)
146        self.plex.metricSetMaximumMagnitude(h_max)
147        self.plex.metricSetMaximumAnisotropy(a_max)
148        self.plex.metricSetTargetComplexity(target)
149        self.plex.metricSetNormalizationOrder(p)
150
151        self.assertFalse(self.plex.metricIsIsotropic())
152        self.assertFalse(self.plex.metricRestrictAnisotropyFirst())
153        self.assertEqual(self.plex.metricGetMinimumMagnitude(), h_min)
154        self.assertEqual(self.plex.metricGetMaximumMagnitude(), h_max)
155        self.assertEqual(self.plex.metricGetMaximumAnisotropy(), a_max)
156        self.assertEqual(self.plex.metricGetTargetComplexity(), target)
157        self.assertEqual(self.plex.metricGetNormalizationOrder(), p)
158
159        metric1 = self.plex.metricCreateUniform(1.0)
160        metric2 = self.plex.metricCreateUniform(2.0)
161        metric = self.plex.metricAverage2(metric1, metric2)
162        metric2.array[:] *= 0.75
163        assert np.allclose(metric.array, metric2.array)
164        metric = self.plex.metricIntersection2(metric1, metric2)
165        assert np.allclose(metric.array, metric1.array)
166        self.plex.metricEnforceSPD(metric)
167        assert np.allclose(metric.array, metric1.array)
168
169    def testAdapt(self):
170        dim = self.plex.getDimension()
171        if dim == 1: return
172        vStart, vEnd = self.plex.getDepthStratum(0)
173        numVertices = vEnd-vStart
174        metric_array = np.zeros([numVertices,dim,dim],dtype=PETSc.ScalarType)
175        for met in metric_array:
176            met[:,:] = np.diag([9]*dim)
177        metric = PETSc.Vec().createWithArray(metric_array)
178        try:
179            newplex = self.plex.adaptMetric(metric,"")
180        except PETSc.Error as exc:
181            if exc.ierr != ERR_SUP: raise
182
183
184# --------------------------------------------------------------------
185
186class BaseTestPlex_2D(BaseTestPlex):
187    DIM = 2
188    CELLS = [[0, 1, 3], [1, 3, 4], [1, 2, 4], [2, 4, 5],
189             [3, 4, 6], [4, 6, 7], [4, 5, 7], [5, 7, 8]]
190    COORDS = [[0.0, 0.0], [0.5, 0.0], [1.0, 0.0],
191              [0.0, 0.5], [0.5, 0.5], [1.0, 0.5],
192              [0.0, 1.0], [0.5, 1.0], [1.0, 1.0]]
193    DOFS = [1, 0, 0]
194
195class BaseTestPlex_3D(BaseTestPlex):
196    DIM = 3
197    CELLS = [[0, 2, 3, 7], [0, 2, 6, 7], [0, 4, 6, 7],
198             [0, 1, 3, 7], [0, 1, 5, 7], [0, 4, 5, 7]]
199    COORDS = [[0., 0., 0.], [1., 0., 0.], [0., 1., 0.], [1., 1., 0.],
200              [0., 0., 1.], [1., 0., 1.], [0., 1., 1.], [1., 1., 1.]]
201    DOFS = [1, 0, 0, 0]
202
203# --------------------------------------------------------------------
204
205class TestPlex_1D(BaseTestPlex, unittest.TestCase):
206    pass
207
208class TestPlex_2D(BaseTestPlex_2D, unittest.TestCase):
209    pass
210
211class TestPlex_3D(BaseTestPlex_3D, unittest.TestCase):
212    pass
213
214class TestPlex_2D_P3(BaseTestPlex_2D, unittest.TestCase):
215    DOFS = [1, 2, 1]
216
217class TestPlex_3D_P3(BaseTestPlex_3D, unittest.TestCase):
218    DOFS = [1, 2, 1, 0]
219
220class TestPlex_3D_P4(BaseTestPlex_3D, unittest.TestCase):
221    DOFS = [1, 3, 3, 1]
222
223class TestPlex_2D_BoxTensor(BaseTestPlex_2D, unittest.TestCase):
224    CELLS = None
225    COORDS = None
226    def setUp(self):
227        self.plex = PETSc.DMPlex().createBoxMesh([3,3], simplex=False)
228
229class TestPlex_3D_BoxTensor(BaseTestPlex_3D, unittest.TestCase):
230    CELLS = None
231    COORDS = None
232    def setUp(self):
233        self.plex = PETSc.DMPlex().createBoxMesh([3,3,3], simplex=False)
234
235try:
236    raise PETSc.Error
237    PETSc.DMPlex().createBoxMesh([2,2], simplex=True, comm=PETSc.COMM_SELF).destroy()
238except PETSc.Error:
239    pass
240else:
241    class TestPlex_2D_Box(BaseTestPlex_2D, unittest.TestCase):
242        CELLS = None
243        COORDS = None
244        def setUp(self):
245            self.plex = PETSc.DMPlex().createBoxMesh([1,1], simplex=True)
246
247    class TestPlex_2D_Boundary(BaseTestPlex_2D, unittest.TestCase):
248        CELLS = None
249        COORDS = None
250        def setUp(self):
251            boundary = PETSc.DMPlex().create(self.COMM)
252            boundary.createSquareBoundary([0., 0.], [1., 1.], [2, 2])
253            boundary.setDimension(self.DIM-1)
254            self.plex = PETSc.DMPlex().generate(boundary)
255
256    class TestPlex_3D_Box(BaseTestPlex_3D, unittest.TestCase):
257        CELLS = None
258        COORDS = None
259        def setUp(self):
260            self.plex = PETSc.DMPlex().createBoxMesh([1,1,1], simplex=True)
261
262    class TestPlex_3D_Boundary(BaseTestPlex_3D, unittest.TestCase):
263        CELLS = None
264        COORDS = None
265        def setUp(self):
266            boundary = PETSc.DMPlex().create(self.COMM)
267            boundary.createCubeBoundary([0., 0., 0.], [1., 1., 1.], [1, 1, 1])
268            boundary.setDimension(self.DIM-1)
269            self.plex = PETSc.DMPlex().generate(boundary)
270
271# --------------------------------------------------------------------
272
273PETSC_DIR = petsc4py.get_config()['PETSC_DIR']
274
275def check_dtype(method):
276    def wrapper(self, *args, **kwargs):
277        if PETSc.ScalarType is PETSc.ComplexType:
278            return
279        else:
280            return method(self, *args, **kwargs)
281    return wrapper
282
283def check_package(method):
284    def wrapper(self, *args, **kwargs):
285        if not PETSc.Sys.hasExternalPackage("hdf5"):
286            return
287        elif self.PARTITIONERTYPE != "simple" and \
288           not PETSc.Sys.hasExternalPackage(self.PARTITIONERTYPE):
289            return
290        else:
291            return method(self, *args, **kwargs)
292    return wrapper
293
294def check_nsize(method):
295    def wrapper(self, *args, **kwargs):
296        if PETSc.COMM_WORLD.size != self.NSIZE:
297            return
298        else:
299            return method(self, *args, **kwargs)
300    return wrapper
301
302class BaseTestPlexHDF5(object):
303    NSIZE = 4
304    NTIMES = 3
305
306    def setUp(self):
307        self.txtvwr = PETSc.Viewer()
308
309    def tearDown(self):
310        if not PETSc.COMM_WORLD.rank:
311            if os.path.exists(self.outfile()):
312                os.remove(self.outfile())
313            if os.path.exists(self.tmp_output_file()):
314                os.remove(self.tmp_output_file())
315        self.txtvwr = None
316
317    def _name(self):
318        return "%s_outformat-%s_%s" % (self.SUFFIX,
319                                       self.OUTFORMAT,
320                                       self.PARTITIONERTYPE)
321
322    def infile(self):
323        return os.path.join(PETSC_DIR, "share/petsc/datafiles/",
324                            "meshes/blockcylinder-50.h5")
325
326    def outfile(self):
327        return os.path.join("./temp_test_dmplex_%s.h5" % self._name())
328
329    def informat(self):
330        return PETSc.Viewer.Format.HDF5_XDMF
331
332    def outformat(self):
333        d = {"hdf5_petsc": PETSc.Viewer.Format.HDF5_PETSC,
334             "hdf5_xdmf": PETSc.Viewer.Format.HDF5_XDMF}
335        return d[self.OUTFORMAT]
336
337    def partitionerType(self):
338        d = {"simple": PETSc.Partitioner.Type.SIMPLE,
339             "ptscotch": PETSc.Partitioner.Type.PTSCOTCH,
340             "parmetis": PETSc.Partitioner.Type.PARMETIS}
341        return d[self.PARTITIONERTYPE]
342
343    def ref_output_file(self):
344        return os.path.join(PETSC_DIR, "src/dm/impls/plex/tutorials/",
345                            "output/ex5_%s.out" % self._name())
346
347    def tmp_output_file(self):
348        return os.path.join("./temp_test_dmplex_%s.out" % self._name())
349
350    def outputText(self, msg, comm):
351        if not comm.rank:
352            with open(self.tmp_output_file(), 'a') as f:
353                f.write(msg)
354
355    def outputPlex(self, plex):
356        self.txtvwr.createASCII(self.tmp_output_file(),
357                                mode='a', comm=plex.comm)
358        plex.view(viewer=self.txtvwr)
359        self.txtvwr.destroy()
360
361    @check_dtype
362    @check_package
363    @check_nsize
364    def testViewLoadCycle(self):
365        grank = PETSc.COMM_WORLD.rank
366        for i in range(self.NTIMES):
367            if i == 0:
368                infname = self.infile()
369                informt = self.informat()
370            else:
371                infname = self.outfile()
372                informt = self.outformat()
373            if self.HETEROGENEOUS:
374                mycolor = (grank > self.NTIMES - i)
375            else:
376                mycolor = 0
377            try:
378                import mpi4py
379            except ImportError:
380                self.skipTest('mpi4py') # throws special exception to signal test skip
381            mpicomm = PETSc.COMM_WORLD.tompi4py()
382            comm = PETSc.Comm(comm=mpicomm.Split(color=mycolor, key=grank))
383            if mycolor == 0:
384                self.outputText("Begin cycle %d\n" % i, comm)
385                plex = PETSc.DMPlex()
386                vwr = PETSc.ViewerHDF5()
387                # Create plex
388                plex.create(comm=comm)
389                plex.setName("DMPlex Object")
390                # Load data from XDMF into dm in parallel
391                vwr.create(infname, mode='r', comm=comm)
392                vwr.pushFormat(format=informt)
393                plex.load(viewer=vwr)
394                plex.setOptionsPrefix("loaded_")
395                plex.setFromOptions()
396                vwr.popFormat()
397                vwr.destroy()
398                self.outputPlex(plex)
399                # Test DM is indeed distributed
400                flg = plex.isDistributed()
401                self.outputText("Loaded mesh distributed? %s\n" %
402                                str(flg).upper(), comm)
403                # Interpolate
404                plex.interpolate()
405                plex.setOptionsPrefix("interpolated_")
406                plex.setFromOptions()
407                self.outputPlex(plex)
408                # Redistribute
409                part = plex.getPartitioner()
410                part.setType(self.partitionerType())
411                _ = plex.distribute(overlap=0)
412                plex.setOptionsPrefix("redistributed_")
413                plex.setFromOptions()
414                self.outputPlex(plex)
415                # Save redistributed dm to XDMF in parallel
416                vwr.create(self.outfile(), mode='w', comm=comm)
417                vwr.pushFormat(format=self.outformat())
418                plex.view(viewer=vwr)
419                vwr.popFormat()
420                vwr.destroy()
421                # Destroy plex
422                plex.destroy()
423                self.outputText("End   cycle %d\n--------\n" % i, comm)
424            PETSc.COMM_WORLD.Barrier()
425        # Check that the output is identical to that of plex/tutorial/ex5.c.
426        self.assertTrue(filecmp.cmp(self.tmp_output_file(),
427                                    self.ref_output_file(), shallow=False),
428                        'Contents of the files not the same.')
429        PETSc.COMM_WORLD.Barrier()
430
431class BaseTestPlexHDF5Homogeneous(BaseTestPlexHDF5):
432    """Test save on N / load on N."""
433    SUFFIX = 0
434    HETEROGENEOUS = False
435
436class BaseTestPlexHDF5Heterogeneous(BaseTestPlexHDF5):
437    """Test save on N / load on M."""
438    SUFFIX = 1
439    HETEROGENEOUS = True
440
441class TestPlexHDF5PETSCSimpleHomogeneous(BaseTestPlexHDF5Homogeneous,
442                                         unittest.TestCase):
443    OUTFORMAT = "hdf5_petsc"
444    PARTITIONERTYPE = "simple"
445
446"""
447Skipping. PTScotch produces different distributions when run
448in a sequence in a single session.
449
450class TestPlexHDF5PETSCPTScotchHomogeneous(BaseTestPlexHDF5Homogeneous,
451                                           unittest.TestCase):
452    OUTFORMAT = "hdf5_petsc"
453    PARTITIONERTYPE = "ptscotch"
454"""
455
456class TestPlexHDF5PETSCParmetisHomogeneous(BaseTestPlexHDF5Homogeneous,
457                                           unittest.TestCase):
458    OUTFORMAT = "hdf5_petsc"
459    PARTITIONERTYPE = "parmetis"
460
461class TestPlexHDF5XDMFSimpleHomogeneous(BaseTestPlexHDF5Homogeneous,
462                                        unittest.TestCase):
463    OUTFORMAT = "hdf5_xdmf"
464    PARTITIONERTYPE = "simple"
465
466"""
467Skipping. PTScotch produces different distributions when run
468in a sequence in a single session.
469
470class TestPlexHDF5XDMFPTScotchHomogeneous(BaseTestPlexHDF5Homogeneous,
471                                          unittest.TestCase):
472    OUTFORMAT = "hdf5_xdmf"
473    PARTITIONERTYPE = "ptscotch"
474"""
475
476class TestPlexHDF5XDMFParmetisHomogeneous(BaseTestPlexHDF5Homogeneous,
477                                          unittest.TestCase):
478    OUTFORMAT = "hdf5_xdmf"
479    PARTITIONERTYPE = "parmetis"
480
481class TestPlexHDF5PETSCSimpleHeterogeneous(BaseTestPlexHDF5Heterogeneous,
482                                           unittest.TestCase):
483    OUTFORMAT = "hdf5_petsc"
484    PARTITIONERTYPE = "simple"
485
486"""
487Skipping. PTScotch produces different distributions when run
488in a sequence in a single session.
489
490class TestPlexHDF5PETSCPTScotchHeterogeneous(BaseTestPlexHDF5Heterogeneous,
491                                             unittest.TestCase):
492    OUTFORMAT = "hdf5_petsc"
493    PARTITIONERTYPE = "ptscotch"
494"""
495
496class TestPlexHDF5PETSCParmetisHeterogeneous(BaseTestPlexHDF5Heterogeneous,
497                                             unittest.TestCase):
498    OUTFORMAT = "hdf5_petsc"
499    PARTITIONERTYPE = "parmetis"
500
501class TestPlexHDF5XDMFSimpleHeterogeneous(BaseTestPlexHDF5Heterogeneous,
502                                          unittest.TestCase):
503    OUTFORMAT = "hdf5_xdmf"
504    PARTITIONERTYPE = "simple"
505
506class TestPlexHDF5XDMFPTScotchHeterogeneous(BaseTestPlexHDF5Heterogeneous,
507                                            unittest.TestCase):
508    OUTFORMAT = "hdf5_xdmf"
509    PARTITIONERTYPE = "ptscotch"
510
511class TestPlexHDF5XDMFParmetisHeterogeneous(BaseTestPlexHDF5Heterogeneous,
512                                            unittest.TestCase):
513    OUTFORMAT = "hdf5_xdmf"
514    PARTITIONERTYPE = "parmetis"
515
516# --------------------------------------------------------------------
517
518if __name__ == '__main__':
519    unittest.main()
520