slepc-main 2024-12-17
Report Typos and Errors

DSGetMat

Returns a sequential dense Mat object containing the requested matrix.

Synopsis

#include "slepcds.h" 
PetscErrorCode DSGetMat(DS ds,DSMatType m,Mat *A)
Not Collective

Input Parameters

ds  - the direct solver context
m  - the requested matrix

Output Parameter

A  - Mat object

Notes

The returned Mat has sizes equal to the current DS dimensions (nxm), and contains the values that would be obtained with DSGetArray() (not DSGetArrayReal()). If the DS was truncated, then the number of rows is equal to the dimension prior to truncation, see DSTruncate(). The communicator is always PETSC_COMM_SELF.

It is implemented with MatDenseGetSubMatrix(), and when no longer needed the user must call DSRestoreMat() which will invoke MatDenseRestoreSubMatrix().

For matrices DS_MAT_T and DS_MAT_D, this function will return a Mat object that cannot be used directly for computations, since it uses compact storage (three and one diagonals for T and D, respectively). In complex scalars, the internal array stores real values, so it is sufficient with 2 columns for T.

See Also

DSRestoreMat(), DSSetDimensions(), DSGetArray(), DSGetArrayReal(), DSTruncate()

Level

advanced

Location

src/sys/classes/ds/interface/dsops.c

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