[petsc-users] Matrix Construction Question

Hong Zhang hzhang at mcs.anl.gov
Thu Jun 30 15:18:48 CDT 2011


Add
 ierr = MatGetOrdering(testMat,MATORDERING_ND,&isrow,&iscol);CHKERRQ(ierr);
before the line
	ierr = MatFactorInfoInitialize(&luinfo);CHKERRQ(ierr);

Somehow, your matrix is numerically singular. With
MATORDERING_NATURAL
I get
[0]PETSC ERROR: Detected zero pivot in LU factorization

even using MATDENSE matrix format, which calls lapack.
With MATORDERING_ND, I get useless inverseMat.

The modified code I used is attached.

Hong

On Thu, Jun 30, 2011 at 2:30 PM, Adam Byrd <adam1.byrd at gmail.com> wrote:
> I'm trying to work through what I need to do, again by practicing with a
> small scale random problem. The general order of events seems to be: create
> a matrix, fill it, assemble it, factor it, then one can use solvers with it.
> When I use MatLUFactor on my matrix before using it with a solver I get this
> error:
> [0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [0]PETSC ERROR: Null argument, when expecting valid pointer!
> [0]PETSC ERROR: Null Object: Parameter # 1!
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: Petsc Release Version 3.1.0, Patch 8, Thu Mar 17 13:37:48
> CDT 2011
> [0]PETSC ERROR: See docs/changes/index.html for recent updates.
> [0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [0]PETSC ERROR: See docs/index.html for manual pages.
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: ./test on a osx-gnu-c named Macintosh-3.local by adambyrd
> Thu Jun 30 15:27:30 2011
> [0]PETSC ERROR: Libraries linked from
> /Users/adambyrd/soft/petsc-3.1-p8/osx-gnu-cpp/lib
> [0]PETSC ERROR: Configure run at Tue Jun 28 12:56:55 2011
> [0]PETSC ERROR: Configure options PETSC_ARCH=osx-gnu-cpp --with-fc=gfortran
> -download-f-blas-lapack=1 --download-mpich=1 --with-scalar-type=complex
> --with-clanguage=c++
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: ISInvertPermutation() line 209 in
> src/vec/is/interface/index.c
> [0]PETSC ERROR: MatLUFactorSymbolic_SeqAIJ() line 306 in
> src/mat/impls/aij/seq/aijfact.c
> [0]PETSC ERROR: MatLUFactorSymbolic() line 2534 in
> src/mat/interface/matrix.c
> [0]PETSC ERROR: MatLUFactor_SeqAIJ() line 945 in
> src/mat/impls/aij/seq/aijfact.c
> [0]PETSC ERROR: MatLUFactor() line 2417 in src/mat/interface/matrix.c
> [0]PETSC ERROR: main() line 62 in WDtest.cpp
> application called MPI_Abort(MPI_COMM_WORLD, 85) - process 0[unset]:
> aborting job:
> application called MPI_Abort(MPI_COMM_WORLD, 85) - process 0
> I don't understand what I'm doing wrong.
> Respectfully,
> Adam
> On Wed, Jun 29, 2011 at 1:26 AM, Matthew Knepley <knepley at gmail.com> wrote:
>>
>> On Tue, Jun 28, 2011 at 10:20 PM, Adam Byrd <adam1.byrd at gmail.com> wrote:
>>>
>>> Matt,
>>>
>>> Alright, that means I need to continue learning how to use
>>> MatSetValues(). With my 6x6 example I tried filling it with four 3x3 sub
>>> matrices, but when I do that I get the error 'sum of local sizes 12 does not
>>> equal global size.' I had 4 processors each calling MatSetValues for their
>>> own 3x3. Graphically, I arranged the nodes 0 1
>>>
>>>                                                                        2 3
>>> where process 0 had global rows 0-2 and global columns 0-2; process 1 had
>>> 0-2, 3-5; process 2 had 3-5, 0-2; and process 3 had 3-5, 3-5. From the
>>> documentation, I think this should be correct, but I'm not sure. Also, which
>>> format would you recommend for storing the matrix?
>>
>> 1) With any error, send the Entire error message.
>> 2) PETSc matrices are divided by rows, not rows and columns, see the
>> manual section. Rows & columns only makes sense for dense matrices
>> 3) You can still set arbitrary blocks no matter how the matrix is divided
>> 4) The error means you tried to set both local and global dimensions, and
>> they do not add up correctly. Just set the global dimensions
>>    Matt
>>
>>>
>>> Jack,
>>>
>>> I'm a summer intern just getting started with this project, so I don't
>>> know all the details yet (I can ask though). I know I need to find the
>>> Green's function which will involve the trace of the inverted Hamiltonian,
>>> as well as the rest of the matrix. I have inquired about avoiding the
>>> inversion altogether, but my instructor doesn't believe there is a way
>>> around it. Once I've worked through the math I want to explore other options
>>> though.
>>>
>>> Respectfully,
>>> Adam
>>>
>>> On Tue, Jun 28, 2011 at 6:08 PM, Matthew Knepley <knepley at gmail.com>
>>> wrote:
>>>>
>>>> On Tue, Jun 28, 2011 at 5:01 PM, Adam Byrd <adam1.byrd at gmail.com> wrote:
>>>>>
>>>>> Actually, it's quite sparse. In the 3600x3600 there are only just 4
>>>>> nonzero entries in each row. This means it's 99.9% empty. My smaller 6x6
>>>>> example is dense, but it's only practice building and manipulating matrices.
>>>>
>>>> Ah, then its easy. Just call MatSetValues() with each block. Then use
>>>> MUMPS to do a sparse direct solve.
>>>>   Matt
>>>>
>>>>>
>>>>> Respectfully,
>>>>> Adam
>>>>>
>>>>> On Tue, Jun 28, 2011 at 5:55 PM, Matthew Knepley <knepley at gmail.com>
>>>>> wrote:
>>>>>>
>>>>>> It sounds like you have a dense matrix (from your example). Is this
>>>>>> true? If so, you should use Elemental (on Google Code).
>>>>>>   Thanks,
>>>>>>      Matt
>>>>>>
>>>>>> On Tue, Jun 28, 2011 at 8:55 AM, Adam Byrd <adam1.byrd at gmail.com>
>>>>>> wrote:
>>>>>>>
>>>>>>> Hi,
>>>>>>> I'm rather new to PETSc and trying to work out the best way to create
>>>>>>> and fill a large sparse matrix distributed over many processors. Currently,
>>>>>>> my goal is to create a 3600x3600 matrix in units of 12x12 blocks with
>>>>>>> several blocks on any given node. I'd like to create the matrix in such a
>>>>>>> way that each node only holds the information in it's handful of blocks and
>>>>>>> not the entire matrix. Eventually, this matrix is to be inverted (I know,
>>>>>>> inversion should be avoided, but as this is a Hamiltonian matrix from which
>>>>>>> I need the Green's function, I'm unaware of a way to forgo carrying out the
>>>>>>> inversion). Additionally, the values will be changed slightly and the matrix
>>>>>>> will be repeatedly inverted. It's structure will remain the same. In order
>>>>>>> to learn how to do this is I am starting with a small 6x6 matrix broken into
>>>>>>> four 3x3 blocks and distributed one block per node. I've been able to create
>>>>>>> a local 3x3 matrix on each node, with it's own values, and with the global
>>>>>>> row/column IDs correctly set to [0, 1, 2] or [3, 4, 5] depending on where
>>>>>>> the block is in the matrix. My problem manifests when I try to create the
>>>>>>> larger matrix from the individual smaller ones. When the matrix is
>>>>>>> constructed I'm trying to use MatSetValues and having each node pass in it's
>>>>>>> 3x3 block. I end up with an error that the sum of local lengths 12x12 does
>>>>>>> not match the global length 6x6. It appears as though this is from passing
>>>>>>> in four 3x3s and the program interpreting that as a 12x12 instead of as a
>>>>>>> 6x6 with the blocks in a grid.
>>>>>>> My question is then: is it possible to fill a matrix as a grid of
>>>>>>> blocks, or can I only fill it in groups of rows or columns? Also, am I
>>>>>>> approaching this problem the correct way, or are there more efficient ways
>>>>>>> of building this matrix with the ultimate goal of inverting it?
>>>>>>> I have included my copy of a modified example if it helps. I do
>>>>>>> apologize if this is answered somewhere in the documentation, I have been
>>>>>>> unable to find a solution.
>>>>>>> Respectfully,
>>>>>>> Adam
>>>>>>
>>>>>>
>>>>>> --
>>>>>> What most experimenters take for granted before they begin their
>>>>>> experiments is infinitely more interesting than any results to which their
>>>>>> experiments lead.
>>>>>> -- Norbert Wiener
>>>>>
>>>>
>>>>
>>>>
>>>> --
>>>> What most experimenters take for granted before they begin their
>>>> experiments is infinitely more interesting than any results to which their
>>>> experiments lead.
>>>> -- Norbert Wiener
>>>
>>
>>
>>
>> --
>> What most experimenters take for granted before they begin their
>> experiments is infinitely more interesting than any results to which their
>> experiments lead.
>> -- Norbert Wiener
>
>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: test.cpp
Type: text/x-c++src
Size: 3138 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20110630/1d6173d1/attachment.cpp>


More information about the petsc-users mailing list