slepc-3.22.1 2024-10-28
Report Typos and Errors

MatCreateTile

Explicitly build a matrix from four blocks, G = [ a*A b*B; c*C d*D ].

Synopsis

#include "slepcsys.h" 
PetscErrorCode MatCreateTile(PetscScalar a,Mat A,PetscScalar b,Mat B,PetscScalar c,Mat C,PetscScalar d,Mat D,Mat *G)
Collective

Input Parameters

A  - matrix for top-left block
a  - scaling factor for block A
B  - matrix for top-right block
b  - scaling factor for block B
C  - matrix for bottom-left block
c  - scaling factor for block C
D  - matrix for bottom-right block
d  - scaling factor for block D

Output Parameter

G  - the resulting matrix

Notes

In the case of a parallel matrix, a permuted version of G is returned. The permutation is a perfect shuffle such that the local parts of A, B, C, D remain in the local part of G for the same process.

Matrix G must be destroyed by the user.

The blocks can be of different type. They can be either ConstantDiagonal, or a standard type such as AIJ, or any other type provided that it supports the MatGetRow operation. The type of the output matrix will be the same as the first block that is not ConstantDiagonal (checked in the A,B,C,D order).

See Also

MatCreateNest()

Level

developer

Location

src/sys/mat/matutil.c

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