[petsc-dev] GAMG error with MKL

Jeff Hammond jeff.science at gmail.com
Tue Jul 3 22:31:53 CDT 2018


On Tue, Jul 3, 2018 at 4:35 PM Mark Adams <mfadams at lbl.gov> wrote:

> On Tue, Jul 3, 2018 at 1:00 PM Richard Tran Mills <rtmills at anl.gov> wrote:
>
>> Hi Mark,
>>
>> I'm glad to see you trying out the AIJMKL stuff. I think you are the
>> first person trying to actually use it, so we are probably going to expose
>> some bugs and also some performance issues. My somewhat limited testing has
>> shown that the MKL sparse routines often perform worse than our own
>> implementations in PETSc.
>>
>
> My users just want OpenMP.
>
>

Why not just add OpenMP to PETSc? I know certain developers hate it, but it
is silly to let a principled objection stand in the way of enabling users
if that would deliver the best performance for NERSC users.

Jeff

We need to systematically look at this and let the MKL team know where
>> performance is lagging.
>>
>> I'm able to run some GAMG examples with AIJMKL matrices if I use type
>> MPIAIJ for my "parallel" matrices but SEQAIJMKL for the "sequential"
>> matrices that MPIAIJ uses. Does your example work if you provide the option
>> "-mat_seqaij_type seqaijmkl"?
>>
>
> Humm, maybe that is a better way to do it. I'll try it. I've been trying
> to keep the MPI matrix as an MKL matrix but I guess that is not necessary.
>
>
>>
>> --Richard
>>
>> On Tue, Jul 3, 2018 at 6:28 AM, Mark Adams <mfadams at lbl.gov> wrote:
>>
>>> GAMG drills into AIJ data structures and will need to be fixed up to
>>> work with MKL matrices, I guess, but it is failing now from a logic error.
>>>
>>> This example works with one processor but fails with 2 (appended). The
>>> code looks like this:
>>>
>>>     ierr =
>>> PetscObjectBaseTypeCompare((PetscObject)Amat,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
>>>     ierr =
>>> PetscObjectBaseTypeCompare((PetscObject)Amat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
>>>     if (isseqaij) {
>>>        ....
>>>     } else if (ismpiaij) {
>>>       .....
>>>       ierr = MatMPIAIJGetSeqAIJ(Amat,&Daij,&Oaij,&garray);CHKERRQ(ierr);
>>>
>>> The failure below is from this call to MatMPIAIJGetSeqAIJ, on this line:
>>>
>>> ierr =
>>> PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
>>>
>>> The difference is PetscObjectTypeCompare vs PetscObject*Base*TypeCompare.
>>> GAMG could use PetscObjectTypeCompare or MatMPIAIJGetSeqAIJ could use
>>> PetscObject*Base*TypeCompare ...
>>>
>>> Mark
>>>
>>> srun -n 2 ./ex19 -pc_type gamg -snes_monitor_short -ksp_monitor_short
>>>
>>>
>>> lid velocity = 0.0625, prandtl # = 1., grashof # = 1.
>>>
>>>
>>>   0 SNES Function norm 0.239155
>>>
>>>
>>> [0]PETSC ERROR: --------------------- Error Message
>>> --------------------------------------------------------------
>>>
>>> [0]PETSC ERROR: No support for this operation for this object type
>>>
>>>
>>> [0]PETSC ERROR: This function requires a MATMPIAIJ matrix as input
>>>
>>>
>>> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html
>>> for trouble shooting.
>>>
>>> [0]PETSC ERROR: Petsc Development GIT revision: v3.9.2-825-g3a11c7608d
>>> GIT Date: 2018-07-01 06:15:09 +0200
>>>
>>> [0]PETSC ERROR:
>>> /global/u2/m/madams/petsc_install/petsc/src/snes/examples/tutorials/./ex19
>>> on a  named nid02516 by madams Tue Jul  3 04:58:13 2018
>>>
>>> [0]PETSC ERROR: Configure options --known-level1-dcache-size=32768
>>> --known-level1-dcache-linesize=64 --known-level1-dcache-assoc=8
>>> --known-sizeof-char=1 --known-sizeof-void-p=8 --known-sizeof-short=2
>>> --known-sizeof-int=4 --known-sizeof-long=8 --known-sizeof-long-long=8
>>> --known-sizeof-float=4 --known-sizeof-double=8 --known-sizeof-size_t=8
>>> --known-bits-per-byte=8 --known-memcmp-ok=1 --known-sizeof-MPI_Comm=4
>>> --known-sizeof-MPI_Fint=4 --known-mpi-long-double=1 --known-mpi-int64_t=1
>>> --known-mpi-c-double-complex=1 --known-mklspblas-supports-zero-based=0
>>> --known-has-attribute-aligned=1 --with-cc=cc --with-cxx=CC --with-fc=ftn
>>> COPTFLAGS="  -g -O1 -mkl -static-intel" CXXOPTFLAGS="-g -O1 -mkl
>>> -static-intel" FOPTFLAGS="  -g -O1 -mkl -static-intel" --download-metis=1
>>> --with-hypre-dir=/global/homes/m/madams/tmp/hypre32-2.14.0
>>> --download-parmetis=1
>>> LIBS="-L/opt/intel/compilers_and_libraries_2018.1.163/linux/mkl/lib/intel64
>>> -lmkl_intel_lp64 -lmkl_intel_thread -lmkl_core -liomp5 -lpthread"
>>> --with-blaslapack-include=/opt/intel/compilers_and_libraries_2018.1.163/linux/mkl/include
>>> --with-debugging=0 --with-mpiexec=srun --with-batch=1
>>> --known-mpi-shared-libraries=1 --known-64-bit-blas-indices=0
>>> --with-64-bit-indices=0 PETSC_ARCH=arch-cori-knl-opt-intel-omp
>>> --with-openmp=1 --download-p4est=0 --with-x=0
>>> --prefix=/global/homes/m/madams/petsc_install/petsc-cori-knl-opt-intel-omp
>>> PETSC_DIR=/global/homes/m/madams/petsc_install/petsc
>>> [0]PETSC ERROR: #1 MatMPIAIJGetSeqAIJ() line 4328 in
>>> /global/u2/m/madams/petsc_install/petsc/src/mat/impls/aij/mpi/mpiaij.c
>>> [0]PETSC ERROR: #2 PCGAMGCreateGraph() line 118 in
>>> /global/u2/m/madams/petsc_install/petsc/src/ksp/pc/impls/gamg/util.c
>>> [0]PETSC ERROR: #3 PCGAMGGraph_AGG() line 832 in
>>> /global/u2/m/madams/petsc_install/petsc/src/ksp/pc/impls/gamg/agg.c
>>> [0]PETSC ERROR: #4 PCSetUp_GAMG() line 517 in
>>> /global/u2/m/madams/petsc_install/petsc/src/ksp/pc/impls/gamg/gamg.c
>>> [0]PETSC ERROR: #5 PCSetUp() line 932 in
>>> /global/u2/m/madams/petsc_install/petsc/src/ksp/pc/interface/precon.c
>>> [0]PETSC ERROR: #6 KSPSetUp() line 381 in
>>> /global/u2/m/madams/petsc_install/petsc/src/ksp/ksp/interface/itfunc.c
>>> [0]PETSC ERROR: #7 KSPSolve() line 612 in
>>> /global/u2/m/madams/petsc_install/petsc/src/ksp/ksp/interface/itfunc.c
>>> [0]PETSC ERROR: #8 SNESSolve_NEWTONLS() line 224 in
>>> /global/u2/m/madams/petsc_install/petsc/src/snes/impls/ls/ls.c
>>> [0]PETSC ERROR: #9 SNESSolve() line 4350 in
>>> /global/u2/m/madams/petsc_install/petsc/src/snes/interface/snes.c
>>> [0]PETSC ERROR: #10 main() line 161 in
>>> /global/homes/m/madams/petsc_install/petsc/src/snes/examples/tutorials/ex19.c
>>> [0]PETSC ERROR: PETSc Option Table entries:
>>> [0]PETSC ERROR: -ksp_monitor_short
>>> [0]PETSC ERROR: -mat_type aijmkl
>>> [0]PETSC ERROR: -options_left
>>> [0]PETSC ERROR: -pc_type gamg
>>> [0]PETSC ERROR: -snes_monitor_short
>>>
>>>
>> --
Jeff Hammond
jeff.science at gmail.com
http://jeffhammond.github.io/
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20180703/f3618c6c/attachment.html>


More information about the petsc-dev mailing list