1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD>
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>MatBDiagGetData</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<A NAME="MatBDiagGetData"><H1>MatBDiagGetData</H1></A>
Gets the data for the block diagonal matrix format. For the parallel case, this returns information for the local submatrix.
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
PetscErrorCode PETSCMAT_DLLEXPORT MatBDiagGetData(Mat mat,PetscInt *nd,PetscInt *bs,PetscInt *diag[],PetscInt *bdlen[],PetscScalar ***diagv)
</PRE>
<H3><FONT COLOR="#CC3333">Input Parameters</FONT></H3>
<DT><B>mat </B> -the matrix, stored in block diagonal format.
<br>
<P>
Not Collective
<P>
<H3><FONT COLOR="#CC3333">Output Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>m </B></TD><TD>- number of rows
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>n </B></TD><TD>- number of columns
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>nd </B></TD><TD>- number of block diagonals
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>bs </B></TD><TD>- each element of a diagonal is an bs x bs dense matrix
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>bdlen </B></TD><TD>- array of total block lengths of block diagonals
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>diag </B></TD><TD>- optional array of block diagonal numbers (length nd).
For a matrix element A[i,j], where i=row and j=column, the
diagonal number is
</TD></TR>
<pre>
diag = i/bs - j/bs (integer division)
</pre>
Set diag=<A HREF="../Sys/PETSC_NULL.html#PETSC_NULL">PETSC_NULL</A> on input for PETSc to dynamically allocate memory as
needed (expensive).
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>diagv </B></TD><TD>- pointer to actual diagonals (in same order as diag array),
</TD></TR></TABLE>
<P>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
See the users manual for further details regarding this storage format.
<P>
<H3><FONT COLOR="#CC3333">Keywords</FONT></H3>
matrix, block, diagonal, get, data
<BR>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../Mat/MatCreateSeqBDiag.html#MatCreateSeqBDiag">MatCreateSeqBDiag</A>(), <A HREF="../Mat/MatCreateMPIBDiag.html#MatCreateMPIBDiag">MatCreateMPIBDiag</A>()
<BR><P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>advanced
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/mat/impls/bdiag/mpi/mpibdiag.c.html#MatBDiagGetData">src/mat/impls/bdiag/mpi/mpibdiag.c</A>
<BR><A HREF="./index.html">Index of all Mat routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>
|