petsc-3.12.5 2020-03-29
Report Typos and Errors

MatTransposeMatMult

Performs Matrix-Matrix Multiplication C=A^T*B.

Synopsis

#include "petscmat.h" 
PetscErrorCode MatTransposeMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
Neighbor-wise Collective on Mat

Input Parameters

A - the left matrix
B - the right matrix
scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use PETSC_DEFAULT if not known

Output Parameters

C -the product matrix

Notes

C will be created if MAT_INITIAL_MATRIX and must be destroyed by the user with MatDestroy().

MAT_REUSE_MATRIX can only be used if the matrices A and B have the same nonzero pattern as in the previous call

To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value actually needed.

This routine is currently implemented for pairs of AIJ matrices and pairs of SeqDense matrices and classes which inherit from SeqAIJ. C will be of same type as the input matrices.

See Also

MatTransposeMatMultSymbolic(), MatTransposeMatMultNumeric(), MatMatMult(), MatMatTransposeMult(), MatPtAP()

Level

intermediate

Location

src/mat/interface/matrix.c

Examples

src/ksp/ksp/examples/tutorials/ex72.c.html

Implementations

MatTransposeMatMult_MPIAIJ_MPIAIJ in src/mat/impls/aij/mpi/mpimatmatmult.c
MatTransposeMatMult_MPIAIJ_MPIDense in src/mat/impls/aij/mpi/mpimattransposematmult.c
MatTransposeMatMult_SeqAIJMKL_SeqAIJMKL_SpMV2 in src/mat/impls/aij/seq/aijmkl/aijmkl.c
MatTransposeMatMult_SeqAIJ_SeqAIJ in src/mat/impls/aij/seq/matmatmult.c
MatTransposeMatMult_SeqAIJ_SeqDense in src/mat/impls/aij/seq/mattransposematmult.c
MatTransposeMatMult_MPIDense_MPIDense in src/mat/impls/dense/mpi/mpidense.c
MatTransposeMatMult_SeqDense_SeqDense in src/mat/impls/dense/seq/dense.c

Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages