ngsPETSc.utils.fenicsx
======================

.. py:module:: ngsPETSc.utils.fenicsx

.. autoapi-nested-parse::

   This module contains all the functions related to wrapping NGSolve meshes to FEniCSx
   We adopt the same docstring conventiona as the FEniCSx project, since this part of
   the package will only be used in combination with FEniCSx.



Attributes
----------

.. autoapisummary::

   ngsPETSc.utils.fenicsx._ngs_to_cells


Classes
-------

.. autoapisummary::

   ngsPETSc.utils.fenicsx.GeometricModel


Functions
---------

.. autoapisummary::

   ngsPETSc.utils.fenicsx._dim_to_element_wrapper
   ngsPETSc.utils.fenicsx._numpy_trim
   ngsPETSc.utils.fenicsx.extract_element_tags


Module Contents
---------------

.. py:data:: _ngs_to_cells

.. py:function:: _dim_to_element_wrapper(ngmesh: Any) -> dict[int, Any]

   Convenience wrapper to extract elements from a NetGen mesh based on topological dimension


.. py:function:: _numpy_trim(array: numpy.typing.ArrayLike, dtype: numpy.typing.DTypeLike) -> numpy.ndarray

   Trim zeros from a two dimensional numpy array along axis 1
   and return it as a specified dtype.

   Args:
       array: The input array to trim.
       dtype: The desired output dtype.
   Returns:
       The trimmed array.


.. py:class:: GeometricModel(geo, comm: mpi4py.MPI.Comm, comm_rank: int = 0)

   This class is used to wrap a Netgen geometric model to a DOLFINx mesh.
   Args:
           geo: The Netgen model
           comm: The MPI communicator to use for mesh creation


   .. py:attribute:: geo


   .. py:attribute:: comm


   .. py:attribute:: comm_rank
      :value: 0



   .. py:method:: model_to_mesh(hmax: float, gdim: int = 2, partitioner: Callable[[mpi4py.MPI.Comm, int, int, dolfinx.cpp.graph.AdjacencyList_int32], dolfinx.cpp.graph.AdjacencyList_int32] | None = None, transform: Any = None, routine: Any = None, meshing_options: dict[str, Any] | None = None) -> tuple[dolfinx.mesh.Mesh, tuple[dolfinx.mesh.MeshTags | None, dolfinx.mesh.MeshTags | None], dict[tuple[int, str], tuple[int, Ellipsis]]]

      Given a NetGen model, take all physical entities of the highest
      topological dimension and create the corresponding linear DOLFINx mesh.

      Args:
          hmax: The maximum diameter of the elements in the triangulation
          model: The NetGen model
          gdim: Geometrical dimension of the mesh
          partitioner: Function that computes the parallel
              distribution of cells across MPI ranks
          transform: PETSc DMPLEX Transformation to be applied to the mesh
          routine: Function to be applied to the mesh after generation
              takes as plan the mesh and the NetGen model and returns the
              same objects after the routine has been applied.

      Returns:
          A DOLFINx mesh for the given NetGen model. It also extracts cell tags,
          facet tags and a mapping from the NetGen label to the corresponding integer marker(s).



   .. py:method:: extract_regions()

      Extract regions from the Netgen mesh.



   .. py:method:: extract_linear_mesh(gdim: int, partitioner: Callable[[mpi4py.MPI.Comm, int, int, dolfinx.cpp.graph.AdjacencyList_int32], dolfinx.cpp.graph.AdjacencyList_int32] | None = None) -> tuple[dolfinx.mesh.MeshTags, dolfinx.mesh.MeshTags]

      Extract a DOLFINx mesh (and correpsonding cell and facet tags) from the Netgen mesh.

      Args:
          gdim: Geometric dimension of the mesh
          partitioner: Function that computes the parallel distribution of cells across MPI ranks

      Note:
          This function updates the `self._mesh` to be in sync with the NetGen model.

      Returns:
          The cell and facet tags of the DOLFINx mesh.



   .. py:method:: curveField(order: int, permutation_tol: float = 1e-08)

      This method returns a curved mesh as a Firedrake function.

      :arg order: the order of the curved mesh.
      :arg permutation_tol: tolerance used to construct the permutation of the reference element.



   .. py:method:: refineMarkedElements(dim: int, elements: numpy.typing.NDArray[numpy.int32], netgen_flags: dict | None = None)

      Refine mesh based on marked elements.



.. py:function:: extract_element_tags(comm_rank: int, ngmesh, dolfinx_mesh: dolfinx.mesh.Mesh, dim: int) -> dolfinx.mesh.MeshTags

   Extract element tags from a Netgen mesh (on a given MPI rank) and distribute them onto
   the corresponding DOLFINx mesh.

   Args:
       comm_rank: The MPI rank to extract the element from.
       ngmesh: The Netgen mesh object.
       dolfinx_mesh: The DOLFINx mesh to which the facet tags will be distributed to.
       dim: The topological dimension of the entities to extract.


