Source code for dune.mmesh._grids

"""
.. module:: _grids
"""

from __future__ import absolute_import, division, print_function, unicode_literals

import sys
import logging
logger = logging.getLogger(__name__)

[docs]def mmesh(constructor, dimgrid=None, **parameters): """Create an MMesh grid. Args: constructor: Grid constructor, e.g. (dune.grid.reader.gmsh, 'grid.msh') or dune.grid.cartesianDomain([0,0],[1,1],[10,10]). dimgrid (int, optional): The world dimension (must be 2 or 3). Might be guessed by constructor. **parameters: Additional parameters. Returns: An MMesh grid instance. Examples: Creating simple structured grid:: from dune.grid import cartesianDomain from dune.mmesh import mmesh grid = mmesh(cartesianDomain([0,0], [1,1], [10,10])) Using a .msh file:: from dune.grid import reader from dune.mmesh import mmesh grid = mmesh((reader.gmsh, "grid.msh"), 2) igrid = grid.hierarchicalGrid.interfaceGrid where grid.msh can be generated by `gmsh <https://pypi.org/project/gmsh/>`_ [GR09]_. For instance:: name = "grid.msh" h = 0.1 hf = 0.01 import gmsh gmsh.initialize() gmsh.option.setNumber("Mesh.MshFileVersion", 2.2) gmsh.model.add("name") geo = gmsh.model.geo p1 = geo.addPoint(0, 0, 0, h, 1) p2 = geo.addPoint(1, 0, 0, h, 2) p3 = geo.addPoint(1, 1, 0, h, 3) p4 = geo.addPoint(0, 1, 0, h, 4) l1 = geo.addLine(p1, p2, 1) l2 = geo.addLine(p2, p3, 2) l3 = geo.addLine(p3, p4, 3) l4 = geo.addLine(p4, p1, 4) p5 = geo.addPoint(0.25, 0.5, 0, hf, 5) p6 = geo.addPoint(0.75, 0.5, 0, hf, 6) lf = geo.addLine(p5, p6, 10) geo.addCurveLoop([l1, l2, l3, l4], 1) geo.addPlaneSurface([1], 0) geo.synchronize() gmsh.model.mesh.embed(1, [lf], 2, 0) gmsh.model.mesh.generate(dim=2) gmsh.write(name) gmsh.finalize() .. [GR09] C. Geuzaine and J.-F. Remacle. Gmsh: a three-dimensional finite element mesh generator with built-in pre- and post-processing facilities. International Journal for Numerical Methods in Engineering 79(11), pp. 1309-1331, 2009. """ if not dimgrid: from dune.grid.grid_generator import getDimgrid dimgrid = getDimgrid(constructor) if not (2 <= dimgrid and dimgrid <= 3): raise KeyError("Parameter error in MMesh with dimgrid=" + str(dimgrid) + ": dimgrid has to be either 2 or 3") if dimgrid == 2: import dune.mmesh._grid2d as gridModule import dune.mmesh._interfacegrid2d as igridModule elif dimgrid == 3: import dune.mmesh._grid3d as gridModule import dune.mmesh._interfacegrid3d as igridModule gridView = gridModule.LeafGrid(gridModule.reader(constructor)) # in case of a cartesian domain store if old or new boundary ids was used # this can be removed in later version - it is only used in dune-fem # to give a warning that the boundary ids for the cartesian domains have changed try: gridView.hierarchicalGrid._cartesianConstructionWithIds = constructor.boundaryWasSet except AttributeError: pass interfaceGrid = lambda self: self.interfaceHierarchicalGrid.leafView gridModule.HierarchicalGrid.interfaceGrid = property(interfaceGrid) return gridView
grid_registry = { "MMesh" : mmesh } if __name__ == "__main__": import doctest doctest.testmod(optionflags=doctest.ELLIPSIS)