slepc-3.15.0 2021-03-31
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
Notes
The Mat is created with sizes equal to the current DS dimensions (nxm),
then it is filled with 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.
When no longer needed, the user can either destroy the matrix or call
DSRestoreMat(). The latter will copy back the modified values.
See Also
DSRestoreMat(), DSSetDimensions(), DSGetArray(), DSGetArrayReal(), DSTruncate()
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