petsc-3.13.6 2020-09-29
Report Typos and Errors

MatGalerkin

Constructs the coarse grid problem via Galerkin projection.

Synopsis

#include "petscmat.h" 
PetscErrorCode  MatGalerkin(Mat restrct, Mat dA, Mat interpolate, MatReuse reuse, PetscReal fill, Mat *A)
If the interpolation and restriction operators are the same, uses MatPtAP. If they are not the same, use MatMatMatMult.

Once the coarse grid problem is constructed, correct for interpolation operators that are not of full rank, which can legitimately happen in the case of non-nested geometric multigrid.

Input Parameters

restrct - restriction operator
dA - fine grid matrix
interpolate - interpolation operator
reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
fill - expected fill, use PETSC_DEFAULT if you do not have a good estimate

Output Parameters

A -the Galerkin coarse matrix

Options Database Key

-pc_mg_galerkin <both,pmat,mat,none> -

See Also

MatPtAP(), MatMatMatMult()

Level

developer

Location

src/mat/interface/matrix.c
Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages