MatMerge_SeqsToMPI

Matthew Knepley knepley at gmail.com
Tue May 20 19:30:06 CDT 2008


On Tue, May 20, 2008 at 7:15 PM, Waad Subber <w_subber at yahoo.com> wrote:
> Thank you Matt and Barry,
>
> The system I am trying to solve is the interface problem in iterative
> substructuring DDM. Where A_i represents [R_i^T*S_i*R_i] and f_i is
> [R_i^T*g_i].
>
> Each process  constructs the local Schur complement matrix (S_i) , the
> restriction matrix(R_i) as SeqAIJ and the RHS vector (g_i) as a sequential
> vector.
>
> Now having the Schur complement matrix for each subdomain, I need to solve
> the interface problem (Sum[R_i^T*S_i*R_i])u=Sum[R_i^T*g_i], .. i=1.. to No.
> of process (subdomains) in parallel.

Barry knows much more than me about substructuring, however:
You could form this matrix in parallel, but I thought that involved significant
communication. I would think that an unassembled form would be better.
To do this you would

  1) Construct the VecScatter from the set of local vectors to the global vector

  2) Create a MatShell and put the VecScatter in the user context

  3) For MatMult(), you would

    a) Scatter the input vector into the local copies

    b) Do each local MatMult()

    c) Scatter the result back to the output vector

Some work, but not that complicated.

   Matt

> For the global vector I construct one MPI vector and use VecGetArray ()  for
> each of the sequential vector then use  VecSetValues () to add the values
> into the global MPI vector. That works fine.
>
> However for the global schur complement matix I try the same idea by
> creating one parallel MPIAIJ matrix and using MatGetArray( ) and
> MatSetValues () in order to add the values to the global matrix.
> MatGetArray( ) gives me only the values without indices, so I don't know how
> to add these valuse to the global MPI matrix.
>
> Thanks agin
>
> Waad
>
> Barry Smith <bsmith at mcs.anl.gov> wrote:
>
> On May 20, 2008, at 3:16 PM, Waad Subber wrote:
>
>> Thank you Matt,
>>
>> Any suggestion to solve the problem I am trying to tackle. I want to
>> solve a linear system:
>>
>> Sum(A_i) u= Sum(f_i) , i=1.... to No. of CPUs.
>>
>> Where A_i is a sparse sequential matrix and f_i is a sequential
>> vector. Each CPU has one matrix and one vector of the same size. Now
>> I want to sum up and solve the system in parallel.
>
> Does each A_i have nonzero entries (mostly) associated with one
> part of the matrix? Or does each process have values
> scattered all around the matrix?
>
> In the former case you should simply create one parallel MPIAIJ
> matrix and call MatSetValues() to put the values
> into it. We don't have any kind of support for the later case, perhaps
> if you describe how the matrix entries come about someone
> would have suggestions on how to proceed.
>
> Barry
>
>>
>>
>> Thanks again
>>
>> Waad
>>
>> Matthew Knepley wrote: On Tue, May 20, 2008 at
>> 2:12 PM, Waad Subber wrote:
>> > Hi,
>> >
>> > I am trying to construct a sparse parallel matrix (MPIAIJ) by
>> adding up
>> > sparse sequential matrices (SeqAIJ) from each CPU. I am using
>> >
>> > MatMerge_SeqsToMPI(MPI_Comm comm,Mat seqmat,PetscInt m,PetscInt
>> n,MatReuse
>> > scall,Mat *mpimat)
>> >
>> > to do that. However, when I compile the code I get the following
>> >
>> > undefined reference to `matmerge_seqstompi_'
>> > collect2: ld returned 1 exit status
>> > make: *** [all] Error 1
>> >
>> > Am I using this function correctly ?
>>
>> These have no Fortran bindings right now.
>>
>> Matt
>>
>> > Thanks
>> >
>> > Waad
>> >
>>
>>
>>
>> --
>> 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