1.2. 2D Orthogonal Grid Generation

Here, we will use the in-house orthogonal mesh generator for a simple Cartesian grid.

To run the code, simply type: jupyter nbconvert --to python --execute mesh_ortho_2D.ipynb.

To convert it to a python file (named mesh_ortho_2D.py), simply type: jupyter nbconvert --to python mesh_ortho_2D.ipynb

To run the python file from the terminal, using N processes, simply type: mpiexec -n <N> python mesh_ortho_2D.py

[ ]:
import os
import sys

sys.path.append("../../..")

from pyopensn.mesh import OrthogonalMeshGenerator, KBAGraphPartitioner, PETScGraphPartitioner
from pyopensn.context import UseColor, Finalize

UseColor(False)

1.2.1. List of Nodes

First, we create a list of nodes that are evenly spaced from -L/2 to +L/2.

[ ]:
length = 2.
n_cells = 10
dx = length / n_cells
nodes = [-length/2 + i * dx for i in range(n_cells + 1)]

1.2.2. Mesh and KBA Partition

We use the OrthogonalMeshGenerator and pass the list of nodes per dimension. Here, we pass 2 times the same list of nodes to create a 2D geometry with square cells. Thus, we create a square domain, of side length L, centered on the origin (0,0).

The domain will be partitioned using the KBAGraphPartitioner into 4 subdomains and the locations of the plane cuts along x and y need to be specified.

[ ]:
meshgen = OrthogonalMeshGenerator(
    node_sets=[nodes, nodes],
    partitioner=KBAGraphPartitioner(
        nx=2,
        ny=2,
        xcuts=[0.0],
        ycuts=[0.0],
    ),
)
grid = meshgen.Execute()

1.2.3. Material IDs

When using the in-house OrthogonalMeshGenerator, no material IDs are assigned. The user needs to assign material IDs to all cells. Here, we have a homogeneous domain, so we assign a material ID with value 0 for each cell in the spatial domain.

[ ]:
grid.SetUniformBlockID(0)

1.2.4. Export the mesh

We export to vtu format. The resulting mesh partition is shown below below1

[ ]:
grid.ExportToPVTU("ortho_2D_KBA")

1.2.5. Mesh (again) and Parmetis partition

Now, we partition the mesh using the Parmetis partitioner, accessed through PETScGraphPartitioner.

[ ]:
meshgen2 = OrthogonalMeshGenerator(
    node_sets=[nodes, nodes],
    partitioner=PETScGraphPartitioner(type='parmetis'),
)
grid2 = meshgen2.Execute()
grid2.SetUniformBlockID(0)

1.2.6. Export the mesh

On such a simple regular mesh, both partitioners are giving the same result. The Parmetis partition is shown below below2

[ ]:
grid2.ExportToPVTU("ortho_2D_Parmetis")

1.2.7. Finalize (for Jupyter Notebook only)

In Python script mode, PyOpenSn automatically handles environment termination. However, this automatic finalization does not occur when running in a Jupyter notebook, so explicit finalization of the environment at the end of the notebook is required. Do not call the finalization in Python script mode, or in console mode.

Note that PyOpenSn’s finalization must be called before MPI’s finalization.

[ ]:
from IPython import get_ipython

def finalize_env():
    Finalize()
    MPI.Finalize()

ipython_instance = get_ipython()
if ipython_instance is not None:
    ipython_instance.events.register("post_execute", finalize_env)