[petsc-users] Setting values of a matrix in a blockwise manner

Tim Steinhoff kandanovian at gmail.com
Wed Jun 10 10:31:16 CDT 2015


One can switch on-the-fly?

Say,

we have a block sparse matrix A initilized via MatSetType(A,MATBAIJ)
and subsequent calls of MatSeqBAIJSetPreallocation and
MatSetValuesBlocked

for use with e.g. the built-in Krylov subspace methods and afterwards
changing type via MatSetType(A,MATAIJ)?

Regards
Tim

2015-06-10 17:24 GMT+02:00 Matthew Knepley <knepley at gmail.com>:
> On Wed, Jun 10, 2015 at 10:20 AM, Tim Steinhoff <kandanovian at gmail.com>
> wrote:
>>
>> I wanted to give the multi frontal LU factorization approach a shot as
>> an example for a direct sparse solver.
>
>
> You can just switch types when you use something other than UMFPACK.
>
>   Thanks,
>
>     Matt
>
>>
>> Regards
>> Tim
>>
>> 2015-06-10 17:00 GMT+02:00 Matthew Knepley <knepley at gmail.com>:
>> > On Wed, Jun 10, 2015 at 9:55 AM, Tim Steinhoff <kandanovian at gmail.com>
>> > wrote:
>> >>
>> >> Thanks for the quick reply!
>> >>
>> >> "Yes, however in order to get improved performance, you need type
>> >> MATBAIJ."
>> >>
>> >> I considered MatSetType(A,MATAIJ); i.e. the non-block type since
>> >> UMFPACK seems to require the seqaij type according to the summary
>> >> page. So do I have to refrain from using the more amiable block-type
>> >> if I want to make use of UMFPACK?
>> >
>> >
>> > Yes, I think so. What do you need in UMFPACK?
>> >
>> >   Thanks,
>> >
>> >     Matt
>> >
>> >>
>> >> 2015-06-10 16:45 GMT+02:00 Matthew Knepley <knepley at gmail.com>:
>> >> > On Wed, Jun 10, 2015 at 9:42 AM, Tim Steinhoff
>> >> > <kandanovian at gmail.com>
>> >> > wrote:
>> >> >>
>> >> >> Hi all
>> >> >>
>> >> >> I want to use Petsc to solve some linear systems via the built-in
>> >> >> Krylov
>> >> >> subspace methods as well as by means of UMFPACK.
>> >> >>
>> >> >> The considered matrix is block sparse with blocks of size 6x6.
>> >> >>
>> >> >> Here is what I came up with after taking a look at some of the
>> >> >> examples
>> >> >>
>> >> >> MPI_Comm comm;
>> >> >> Mat A;
>> >> >> PetscInt n = 10000; /* dimension of matrix */
>> >> >> comm = PETSC_COMM_SELF;
>> >> >> MatCreate(comm,&A);
>> >> >> MatSetSizes(Amat,n,n,n,n);
>> >> >> MatSetBlockSize(A,6);
>> >> >> MatSetType(A,MATAIJ); /* UMFPACK compatible format due to comm =
>> >> >> PETSC_COMM_SELF */
>> >> >>
>> >> >> Questions:
>> >> >> 1.
>> >> >> I work on a single node with 2-8 cores. Hence, comm =
>> >> >> PETSC_COMM_SELF;
>> >> >> I
>> >> >> guess. Is it correct in this contect to set
>> >> >> MatSetSizes(Amat,n,n,n,n);
>> >> >> with
>> >> >> 4-times n?
>> >> >
>> >> >
>> >> > Yes.
>> >> >
>> >> >>
>> >> >> 2.
>> >> >> After the above sequence of commands do I have to use something like
>> >> >>   MatSeqAIJSetPreallocation(A,0,d_nnz); /* d_nnz <-> number of
>> >> >> nonzeros
>> >> >> per row */
>> >> >> or is it possible to use
>> >> >>   MatSeqBAIJSetPreallocation(A,6,0,db_nnz); /* db_nnz <-> number of
>> >> >> block
>> >> >> nonzeros per block row */
>> >> >
>> >> >
>> >> > You should use this if using MATBAIJ.
>> >> >
>> >> >>
>> >> >> In any case, is something like
>> >> >>   MatSetValuesBlocked(A,1,idx_r,1,idx_c,myblockvals,INSERT_VALUES);
>> >> >> to fill values of one block into the matrix A ok?
>> >> >
>> >> >
>> >> > Yes, however in order to get improved performance, you need type
>> >> > MATBAIJ.
>> >> >
>> >> >   Thanks,
>> >> >
>> >> >     Matt
>> >> >
>> >> >>
>> >> >> Regards
>> >> >> Tim
>> >> >
>> >> >
>> >> >
>> >> >
>> >> > --
>> >> > 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


More information about the petsc-users mailing list