ngsPETSc.pc
===========

.. py:module:: ngsPETSc.pc

.. autoapi-nested-parse::

   This module contains all the function and class needed to wrap a PETSc Preconditioner in NGSolve



Classes
-------

.. autoapisummary::

   ngsPETSc.pc.PETScPreconditioner
   ngsPETSc.pc.ASMPreconditioner


Functions
---------

.. autoapisummary::

   ngsPETSc.pc.createPETScPreconditioner


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

.. py:class:: PETScPreconditioner(mat, freeDofs, solverParameters=None, optionsPrefix='', nullspace=None, matType='aij')

   Bases: :py:obj:`ngsolve.BaseMatrix`


   This class creates a Netgen/NGSolve BaseMatrix corresponding to a PETSc PC  

   :arg mat: NGSolve Matrix one would like to build the PETSc preconditioner for.

   :arg freeDofs: not constrained degrees of freedom of the finite element space over
   which the BilinearForm corresponding to the matrix is defined.

   :arg solverParameters: parameters to be passed to the KSP solver

   :arg optionsPrefix: special solver options prefix for this specific Krylov solver

   :arg matType: type of sparse matrix, i.e. PETSc sparse: aij, 
   MKL sparse: mklaij or CUDA: aijcusparse



   .. py:attribute:: ngsMat


   .. py:attribute:: vecMap


   .. py:attribute:: lgmap


   .. py:attribute:: petscPreconditioner


   .. py:method:: Shape()

      Shape of the BaseMatrix




   .. py:method:: CreateVector(col)

      Create vector corresponding to the matrix

      :arg col: True if one want a column vector




   .. py:method:: Mult(x, y)

      BaseMatrix multiplication Ax = y
      :arg x: vector we are multiplying
      :arg y: vector we are storeing the result in




   .. py:method:: MultTrans(x, y)

      BaseMatrix multiplication A^T x = y
      :arg x: vector we are multiplying
      :arg y: vector we are storeing the result in




   .. py:method:: setActingDofs(dofs)

      Set the acting dofs of the preconditioner
      :arg dofs: dofs that the preconditioner is acting on



.. py:function:: createPETScPreconditioner(mat, freeDofs, solverParameters)

   Create PETSc PC that can be accessed by NGSolve.
   :arg mat: NGSolve Matrix one would like to build the PETSc preconditioner for.

   :arg freeDofs: not constrained degrees of freedom of the finite element space over
   which the BilinearForm corresponding to the matrix is defined.

   :arg solverParameters: parameters to be passed to the KSP solver



.. py:class:: ASMPreconditioner(mat, freeDofs, solverParameters=None, optionsPrefix='', nullspace=None, matType='aij', blocks=None)

   Bases: :py:obj:`PETScPreconditioner`


   This class creates a Netgen/NGSolve BaseMatrix corresponding to a PETSc ASM PC


   .. py:attribute:: asmpc
      :value: None



   .. py:method:: Mult(x, y)

      BaseMatrix multiplication Ax = y
      :arg x: vector we are multiplying
      :arg y: vector we are storeing the result in




   .. py:method:: MultTrans(x, y)

      BaseMatrix multiplication A^T x = y
      :arg x: vector we are multiplying
      :arg y: vector we are storeing the result in




