[petsc-users] Various Questions Regarding PETSC

Dave May dave.mayhem23 at gmail.com
Wed Jul 17 10:24:26 CDT 2019


Can you please send your code?

It would greatly help Matt (and others looking at this thread)
to understand the performance you see and to offer advice on how (if
possible) to improve it.


On Wed, 17 Jul 2019 at 16:00, Matthew Knepley via petsc-users <
petsc-users at mcs.anl.gov> wrote:

> On Wed, Jul 17, 2019 at 8:51 AM Mohammed Mostafa <
> mo7ammedmostafa at gmail.com> wrote:
>
>> Sorry for the confusion
>> First I fully acknowledge that setting Matrix non-zeros or copying in
>> general is not cheap and memory access pattern can play an important role.
>> So to establish a baseline to compare with, I tried setting the same
>> matrix but in an Eigen Sparse Matrix  and the timings are as follows
>> FillPetscMat_with_MatSetValues              100 1.0 3.8820e+00 1.1
>> 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 23  0  0  0  0 * 96 * 0  0  0  0
>>   0
>> FilEigenMat                                               100 1.0
>> 2.8727e+00 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 18  0  0  0  0  88  0
>>  0  0  0     0
>>
>
> Great. This helps. Two things would help me narrow down what is happening.
>
>   1) Are you timing only the insertion of values, or computation and
> insertion?
>
>   2) Can you tell me how many values are inserted?
>
>   Thanks,
>
>     Matt
>
>
>> I used the same code but simply filled a different Matrix something like
>>
>> for ( i =0; i < nRows;i++)
>> {
>> //
>> .......
>> // Some code to get  j_index, coefvalues
>> // Method1
>> MatSetValues(A, 1, &cell_global_index, nnz_per_row, j_index, coefvalues,
>> INSERT_VALUES);
>>
>> //Method2
>> for ( int k = 0;k < nnz_per_row; k++)
>>      EigenMat.coeffRef(i, j_index[k] ) = coefvalues[k];
>>
>>  }
>> Please note that only one of the two methods is being used at a time.
>> Also, I separately time the code section used to <  j_index,   coefvalues>
>> but simpling disabling both Method1 and Method2.
>> I found the cost to be trivial in comparison to when either one of the
>> methods is used.
>> I used Eigen out of convenience since I used for some vector and tensor
>> arithmetics somewhere else in the code and it may not be the best choice.
>> Since in PetscMatrix we technically fill two matrices: diagonal and
>> off-diagonal so I expected some difference but is that normal or am I
>> missing something. ?
>> Maybe some setting or MatOption I should be using so far this what I have
>> been using
>>
>> MatCreateMPIAIJWithArrays(PETSC_COMM_WORLD, local_size, local_size,
>> PETSC_DETERMINE,
>> PETSC_DETERMINE, ptr, j , v, A);
>>  MatSetOption(A,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
>> MatSetOption(A,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);
>> MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
>> MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
>> MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
>> MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
>> MatSetUp(A);
>> MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
>> MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
>>
>> Thanks,
>> Kamra
>>
>
>
> --
> 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
>
> https://www.cse.buffalo.edu/~knepley/
> <http://www.cse.buffalo.edu/~knepley/>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20190717/58a43aba/attachment.html>


More information about the petsc-users mailing list