[petsc-dev] MatMatMult

Pierre Jolivet pierre.jolivet at enseeiht.fr
Wed Jun 7 01:29:14 CDT 2017


> On 6 Jun 2017, at 11:21 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
> 
>> 
>> On Jun 6, 2017, at 1:52 PM, Pierre Jolivet <Pierre.Jolivet at enseeiht.fr> wrote:
>> 
>> On Tue, 6 Jun 2017 13:44:06 -0500, Hong wrote:
>>> Pierre :
>>> 
>>>> Yes, of course I defined (*C)->ops->matmultnumeric =
>>>> MatMatMultNumeric_MPIBAIJ_MPIDense in
>>>> MatMatMultSymbolic_MPIBAIJ_MPIDense.
>>>> However, the routine MatMatMultSymbolic_MPIBAIJ_MPIDense is never
>>>> reached when calling MatMatMult with scall == MAT_REUSE_MATRIX
>>>> 
>>> (http://www.mcs.anl.gov/petsc/petsc-current/src/mat/interface/matrix.c.html#line9487
>>>> [1], just to be sure I added a dummy printf in
>>>> MatMatMultSymbolic_MPIBAIJ_MPIDense and nothing is displayed with
>>>> MAT_REUSE_MATRIX, MAT_INITIAL_MATRIX works as intended)
>>> 
>>> MatMatMultSymbolic_xxx() is called only for MAT_INITIAL_MATRIX,
>>> during which, it defines
>>> (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIBAIJ_MPIDense;
>> 
>> Sorry, I was not thorough enough. What I meant to say is that I never call MatMatMult with MAT_INITIAL_MATRIX, so MatMatMultSymbolic_xxx is never called.
> 
>   It is not designed to work that way. The reuse is specifically to reuse a matrix obtained with a call to initial_matrix, not to reuse a random matrix the user provides. The reason is that the symbolic stage may/does build additional data structures that are used for the multiples that the user cannot see or build directly.

I’m doing:
  MatCreateDense(PETSC_COMM_WORLD, m, n, M, N, in, &B);
  MatCreateDense(PETSC_COMM_WORLD, m, n, M, N, out, &C);
  MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C);

This works perfectly fine with A of type MPIAIJ, but once again, everything is hardwired for this only type.

According to your answer, Barry, I think the manual should be adjusted accordingly because to me what I’m trying to achieve should be possible as is (the fact that it is possible with MPIAIJ is also puzzling):
"In the special case where matrix B (and hence C) are dense you can create the correctly sized matrix C yourself and then call this routine with MAT_REUSE_MATRIX, rather than first having MatMatMult() create it for you. […]” (http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatMatMult.html <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatMatMult.html>)

>> Is there another way?
>> I cannot use MAT_INITIAL_MATRIX because I want to manage the underlying memory, and I'm guessing MAT_INITIAL_MATRIX will delete whatever is in C and reallocate memory on its own.
> 
>  The memory in C is just a dense array so presumably you want to use some "application" memory for this space? 
> 
>   The simplest way for you to achieve this is to call the MatMatMultSymbolic() then dig out the underly array of the matrix, call PetscFree() on it and put in your own array. Yes, this is slightly hacky, but it will work fine within the current PETSc model. 

Yeah, alright, I already need to tinker with the MatDensePlaceArray/MatDenseResetArray so I can figure something out as well for this.
Thanks,
Pierre

>  Barry
> 
>> 
>> Thanks in advance for your help,
>> Pierre
>> 
>>> Then MatMatMult(A,C,MAT_REUSE_MATRIX,..) calls 
>>> (*(*C)->ops->matmultnumeric)(A,B,*C); (line 9432 in matrix.c)
>>> 
>>> which should go to MatMatMultNumeric_MPIBAIJ_MPIDense.
>>> 
>>> You may follow a debugging process
>>> using petsc/src/mat/examples/tests/ex109.c
>>> 
>>> Are you working on a branch of petsc? If so, I may take a look and
>>> see what is the problem.
>>> 
>>> Hong
>>> 
>>>>>> 2) I'm having trouble when scall == MAT_REUSE_MATRIX. Here,
>>>>> 
>>>>> 
>>>>> 
>>>> 
>>> http://www.mcs.anl.gov/petsc/petsc-current/src/mat/impls/dense/mpi/mpidense.c.html#line1208
>>>>> [2]
>>>>> 
>>>>>> [2] it looks that the numeric part of the MatMatMult (which is
>>>>>> called when scall == MAT_REUSE_MATRIX) is hardwired to this
>>>>>> routine
>>>>> 
>>>>> 
>>>>> 
>>>> 
>>> http://www.mcs.anl.gov/petsc/petsc-current/src/mat/impls/aij/mpi/mpimatmatmult.c.html#line376
>>>>> [3]
>>>>> [3]. 
>>>>> br>
>>>>> 
>>>> 
>>> s.anl.gov/petsc/petsc-current/src/mat/impls/aij/mpi/mpimatmatmult.c.html#line556"
>>>>> rel="noreferrer"
>>>>> 
>>>> 
>>> target="_blank">http://www.mcs.anl.gov/petsc/petsc-current/src/mat/impls/aij/mpi/mpimatmatmult.c.html#line556
>>>>> [2]
>>>>> 
>>>>> 
>>>>> 
>>>> 
>>> http://www.mcs.anl.gov/petsc/petsc-current/src/mat/impls/dense/mpi/mpidense.c.html#line1208
>>>>> [4]
>>>>> [3]
>>>>> 
>>>>> 
>>>>> 
>>>> 
>>> http://www.mcs.anl.gov/petsc/petsc-current/src/mat/impls/aij/mpi/mpimatmatmult.c.html#line376
>>>>> [5]
>>> 
>>> 
>>> 
>>> Links:
>>> ------
>>> [1]
>>> http://www.mcs.anl.gov/petsc/petsc-current/src/mat/interface/matrix.c.html#line9487
>>> [2]
>>> http://www.mcs.anl.gov/petsc/petsc-current/src/mat/impls/dense/mpi/mpidense.c.html#line1208
>>> [3]
>>> http://www.mcs.anl.gov/petsc/petsc-current/src/mat/impls/aij/mpi/mpimatmatmult.c.html#line376
>>> [4]
>>> http://www.mcs.anl.gov/petsc/petsc-current/src/mat/impls/dense/mpi/mpidense.c.html#line1208
>>> [5]
>>> http://www.mcs.anl.gov/petsc/petsc-current/src/mat/impls/aij/mpi/mpimatmatmult.c.html#line376

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20170607/51507d17/attachment.html>


More information about the petsc-dev mailing list