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 61 62 63 64 65 66 67 68 69 70 71
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD>
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>MatCreateSeqBDiag</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<A NAME="MatCreateSeqBDiag"><H1>MatCreateSeqBDiag</H1></A>
Creates a sequential block diagonal matrix.
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
int MatCreateSeqBDiag(MPI_Comm comm,int m,int n,int nd,int bs,const int diag[],PetscScalar *diagv[],Mat *A)
</PRE>
Collective on <A HREF="../Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</A>
<P>
<H3><FONT COLOR="#CC3333">Input Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B><A HREF="../Sys/comm.html#comm">comm</A> </B></TD><TD>- MPI communicator, set to <A HREF="../Sys/PETSC_COMM_SELF.html#PETSC_COMM_SELF">PETSC_COMM_SELF</A>
</TD></TR>
<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 (optional)
</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>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),
if allocated by user. Otherwise, set diagv=<A HREF="../Sys/PETSC_NULL.html#PETSC_NULL">PETSC_NULL</A> on input for PETSc
to control memory allocation.
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Output Parameters</FONT></H3>
<DT><B>A </B> -the matrix
<br>
<P>
<H3><FONT COLOR="#CC3333">Options Database Keys</FONT></H3>
<DT><B>-mat_block_size <bs> </B> -Sets blocksize
<br>
<DT><B>-mat_bdiag_diags <s1,s2,s3,...> </B> -Sets diagonal numbers
<br>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
See the users manual for further details regarding this storage format.
<P>
<H3><FONT COLOR="#CC3333">Fortran Note</FONT></H3>
Fortran programmers cannot set diagv; this value is ignored.
<P>
<P>
<H3><FONT COLOR="#CC3333">Keywords</FONT></H3>
matrix, block, diagonal, sparse
<BR>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../Mat/MatCreate.html#MatCreate">MatCreate</A>(), <A HREF="../Mat/MatCreateMPIBDiag.html#MatCreateMPIBDiag">MatCreateMPIBDiag</A>(), <A HREF="../Mat/MatSetValues.html#MatSetValues">MatSetValues</A>()
<BR><P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>intermediate
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/mat/impls/bdiag/seq/bdiag.c.html#MatCreateSeqBDiag">src/mat/impls/bdiag/seq/bdiag.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>
|