[petsc-dev] building with MKL

Jeff Hammond jeff.science at gmail.com
Sun Jul 1 12:30:29 CDT 2018


I discourage multi-arch builds in the context of KNL. Using Haswell as the
baseline isn’t ideal for KNL uarch.

FLAGS="-xCORE-AVX2 -axMIC-AVX512,COMMON-AVX512 -O2 -g"
CXXOPTFLAGS="-xCORE-AVX2 -axMIC-AVX512,COMMON-AVX512 -O2 -g”

If you really want to do this, then replace COMMON with CORE to specialize
for SKX. There’s not point to using COMMON if you’ve got the MIC path
already.

The other issue is the SKX really should use -xCORE-AVX512
-qopt-zmm-usage=low, and I don’t know how that interoperable with the -ax
options.

Jeff

On Sun, Jul 1, 2018 at 9:49 AM Victor Eijkhout <eijkhout at tacc.utexas.edu>
wrote:

> My petsc installation work with
>
>
> --with-blaslapack-dir=/opt/intel/compilers_and_libraries_2017.4.196/linux/mkl
>
> FWIW
>
> and that seems to be the only thing, but just in case:
>
> Starting configure run at Mon, 30 Apr 2018 15:05:30 +0000
> Configure Options: --configModules=PETSc.Configure
> --optionsModule=config.compilerOptions --with-x=0 -with-pic
> --with-make-np=12 --with-packages-dir=/admin/build/admin/rpms/s\
> tampede2/SOURCES/petsc-packages/externalpackages-3.9
> --with-external-packages-dir=/admin/build/admin/rpms/stampede2/SOURCES/petsc-packages/externalpackages-3.9
> --with-mpi-com\
> pilers=1
> --with-mpi-include=/opt/intel/compilers_and_libraries_2017.4.196/linux/mpi/intel64/include
> --with-mpi-lib=/opt/intel/compilers_and_libraries_2017.4.196/linux/mpi/int\
> el64/lib/release_mt/libmpi.so --with-scalar-type=real
> --with-shared-libraries=1 --with-precision=double --with-chaco=1
> --download-chaco --with-hypre=1 --download-hypre --with\
> -parmetis=1 --download-parmetis --with-metis=1 --download-metis
> --with-plapack=1 --download-plapack --with-spai=1 --download-spai
> --with-sundials=1 --download-sundials --with\
> -hdf5=1 --with-hdf5-dir=/opt/apps/intel17/impi17_0/phdf5/1.8.16/x86_64
> --with-mumps=1 --download-mumps --with-parmetis=1 --download-parmetis
> --with-metis=1 --download-metis -\
> -with-scalapack=1 --download-scalapack --with-blacs=1 --download-blacs
> --with-spooles=1 --download-spooles --with-suitesparse=1
> --download-suitesparse --with-superlu_dist=1 -\
> -download-superlu_dist --with-superlu=1 --download-superlu
> --with-parmetis=1 --download-parmetis --with-metis=1 --download-metis
> --with-zoltan=1 --download-zoltan=1 --downloa\
> d-ptscotch=1 --with-debugging=no
> --with-blaslapack-dir=/opt/intel/compilers_and_libraries_2017.4.196/linux/mkl
> COPTFLAGS="-xCORE-AVX2 -axMIC-AVX512,COMMON-AVX512 -O2 -g" FOPT\
> FLAGS="-xCORE-AVX2 -axMIC-AVX512,COMMON-AVX512 -O2 -g"
> CXXOPTFLAGS="-xCORE-AVX2 -axMIC-AVX512,COMMON-AVX512 -O2 -g”
>
> Victor.
>
>
> On Jul 1, 2018, at 9:17 AM, Mark Adams <mfadams at lbl.gov> wrote:
>
>
>
> On Sun, Jul 1, 2018 at 4:08 AM Karl Rupp <rupp at iue.tuwien.ac.at> wrote:
>
>> Hi Mark,
>>
>> have a look at config/examples/arch-linux-knl.py, which contains on line
>> 20:
>>   '--with-blaslapack-dir='+os.environ['MKLROOT'],
>>
>> It's important that you specify the BLAS library *and* the MKL include
>> directory (either via --with-blaslapack-dir or via a pair of
>> --with-blaslapack-include and --with-blaslapack-lib),
>
>
> I think what I sent was this second option. I nuked, reconfigured and
> tested again. Same thing:
>
> > [0]PETSC ERROR: --------------------- Error Message
> --------------------------------------------------------------
> > [0]PETSC ERROR: Unknown type. Check for miss-spelling or missing package:
>  http://www.mcs.anl.gov/petsc/documentation/installation.html#external
> > [0]PETSC ERROR: Unknown Mat type given: aijmkl
> > [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for
> trouble shooting.
> > [0]PETSC ERROR: Petsc Release Version 3.9.2, unknown
> > [0]PETSC ERROR:
> /global/u2/m/madams/petsc_install/petsc/src/snes/examples/tutorials/./ex19
> on a  named nid02516 by madams Sun Jul  1 07:10:33 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-has-attribute-aligned=1 --with-cc=cc
> --with-cxx=CC --with-fc=ftn COPTFLAGS="  -g -O0 -hcpu=mic-knl
> -qopenmp-simd" CXXOPTFLAGS="-g -O0 -hcpu=mic-knl -qopenmp-simd"
> FOPTFLAGS="  -g -O0 -hcpu=mic-knl -qopenmp-simd" --download-metis=1
> --with-hypre-dir=/global/homes/m/madams/tmp/hypre-2.14.0
> --download-parmetis=1
> --with-blaslapack-lib=/opt/intel/compilers_and_libraries_2018.1.163/linux/mkl/lib/intel64/libmkl_intel_thread.a
> --with-blaslapack-include=/opt/intel/compilers_and_libraries_2018.1.163/linux/mkl/include
> --with-debugging=1 --with-mpiexec=srun --with-batch=1
> --known-mpi-shared-libraries=1 --known-64-bit-blas-indices=0
> --with-64-bit-indices=1 PETSC_ARCH=arch-cori-knl-dbg64-intel-omp
> --with-openmp=1 --download-p4est=0 --with-x=0
> --prefix=/global/homes/m/madams/petsc_install/petsc-cori-knl-dbg64-intel-omp
> PETSC_DIR=/global/homes/m/madams/petsc_install/petsc
> > [0]PETSC ERROR: #1 MatSetType() line 61 in
> /global/u2/m/madams/petsc_install/petsc/src/mat/interface/matreg.c
> > [0]PETSC ERROR: #2 MatSetFromOptions() line 229 in
> /global/u2/m/madams/petsc_install/petsc/src/mat/utils/gcreate.c
> > [0]PETSC ERROR: #3 DMCreateMatrix_DA() line 793 in
> /global/u2/m/madams/petsc_install/petsc/src/dm/impls/da/fdda.c
> > [0]PETSC ERROR: #4 DMCreateMatrix() line 1262 in
> /global/u2/m/madams/petsc_install/petsc/src/dm/interface/dm.c
> > [0]PETSC ERROR: #5 SNESSetUpMatrices() line 646 in
> /global/u2/m/madams/petsc_install/petsc/src/snes/interface/snes.c
> > [0]PETSC ERROR: #6 SNESSetUp_NEWTONLS() line 296 in
> /global/u2/m/madams/petsc_install/petsc/src/snes/impls/ls/ls.c
> > [0]PETSC ERROR: #7 SNESSetUp() line 2908 in
> /global/u2/m/madams/petsc_install/petsc/src/snes/interface/snes.c
> > [0]PETSC ERROR: #8 SNESSolve() line 4300 in
> /global/u2/m/madams/petsc_install/petsc/src/snes/interface/snes.c
> > [0]PETSC ERROR: #9 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: -da_refine 3
> > [0]PETSC ERROR: -ksp_monitor
> > [0]PETSC ERROR: -mat_type aijmkl
> > [0]PETSC ERROR: -options_left
> > [0]PETSC ERROR: -pc_type gamg
> > [0]PETSC ERROR: -snes_monitor_short
> > [0]PETSC ERROR: -snes_view
> > [0]PETSC ERROR: ----------------End of Error Message -------send entire
> error message to petsc-maint at mcs.anl.gov----------
> > Rank 0 [Sun Jul  1 07:10:33 2018] [c1-1c0s5n0] application called
> MPI_Abort(MPI_COMM_WORLD, 86) - process 0
> > srun: error: nid02516: task 0: Aborted
> > srun: Terminating job step 13487457.2
>
>
>
>> otherwise it's not
>> possible to compile the aijmkl code.
>>
>> Best regards,
>> Karli
>>
>>
>>
>>
>> On 06/30/2018 09:55 PM, Mark Adams wrote:
>> > It builds and runs but looks like PETSc does not register aijmkl
>> matrices.
>> >
>> >
>> > ---
>> >  > [0]PETSC ERROR: --------------------- Error Message
>> > --------------------------------------------------------------
>> >  > [0]PETSC ERROR: Unknown type. Check for miss-spelling or missing
>> > package:
>> > http://www.mcs.anl.gov/petsc/documentation/installation.html#external
>> >  > [0]PETSC ERROR: Unknown Mat type given: aijmkl
>> >  > [0]PETSC ERROR: See
>> > http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble
>> shooting.
>> >  > [0]PETSC ERROR: Petsc Release Version 3.9.2, unknown
>> >  > [0]PETSC ERROR:
>> >
>> /global/u2/m/madams/petsc_install/petsc/src/snes/examples/tutorials/./ex19
>>
>> > on a  named nid02516 by madams Sat Jun 30 12:48:10 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-has-attribute-aligned=1 --with-cc=cc --with-cxx=CC --with-fc=ftn
>>
>> > COPTFLAGS="  -g -O0 -hcpu=mic-knl -qopenmp-simd" CXXOPTFLAGS="-g -O0
>> > -hcpu=mic-knl -qopenmp-simd" FOPTFLAGS="  -g -O0 -hcpu=mic-knl
>> > -qopenmp-simd" --download-metis=1
>> > --with-hypre-dir=/global/homes/m/madams/tmp/hypre-2.14.0
>> > --download-parmetis=1
>> >
>> --with-blaslapack-lib=/opt/intel/compilers_and_libraries_2018.1.163/linux/mkl/lib/intel64/libmkl_intel_thread.a
>>
>> > --with-debugging=1 --with-mpiexec=srun --with-batch=1
>> > --known-mpi-shared-libraries=1 --known-64-bit-blas-indices=0
>> > --with-64-bit-indices=1 PETSC_ARCH=arch-cori-knl-dbg64-intel-omp
>> > --with-openmp=1 --download-p4est=0 --with-x=0
>> >
>> --prefix=/global/homes/m/madams/petsc_install/petsc-cori-knl-dbg64-intel-omp
>>
>> > PETSC_DIR=/global/homes/m/madams/petsc_install/petsc
>> >  > [0]PETSC ERROR: #1 MatSetType() line 61 in
>> > /global/u2/m/madams/petsc_install/petsc/src/mat/interface/matreg.c
>> >  > [0]PETSC ERROR: #2 MatSetFromOptions() line 229 in
>> > /global/u2/m/madams/petsc_install/petsc/src/mat/utils/gcreate.c
>> >  > [0]PETSC ERROR: #3 DMCreateMatrix_DA() line 793 in
>> > /global/u2/m/madams/petsc_install/petsc/src/dm/impls/da/fdda.c
>> >  > [0]PETSC ERROR: #4 DMCreateMatrix() line 1262 in
>> > /global/u2/m/madams/petsc_install/petsc/src/dm/interface/dm.c
>> >  > [0]PETSC ERROR: #5 SNESSetUpMatrices() line 646 in
>> > /global/u2/m/madams/petsc_install/petsc/src/snes/interface/snes.c
>> >  > [0]PETSC ERROR: #6 SNESSetUp_NEWTONLS() line 296 in
>> > /global/u2/m/madams/petsc_install/petsc/src/snes/impls/ls/ls.c
>> >  > [0]PETSC ERROR: #7 SNESSetUp() line 2908 in
>> > /global/u2/m/madams/petsc_install/petsc/src/snes/interface/snes.c
>> >  > [0]PETSC ERROR: #8 SNESSolve() line 4300 in
>> > /global/u2/m/madams/petsc_install/petsc/src/snes/interface/snes.c
>> >  > [0]PETSC ERROR: #9 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: -da_refine 3
>> >  > [0]PETSC ERROR: -ksp_monitor
>> >  > [0]PETSC ERROR: -mat_type aijmkl
>> >  > [0]PETSC ERROR: -options_left
>> >  > [0]PETSC ERROR: -pc_type gamg
>> >  > [0]PETSC ERROR: -snes_monitor_short
>> >  > [0]PETSC ERROR: -snes_view
>> >  > [0]PETSC ERROR: ----------------End of Error Message -------se
>> >
>> > On Sat, Jun 30, 2018 at 3:08 PM Mark Adams <mfadams at lbl.gov
>> > <mailto:mfadams at lbl.gov>> wrote:
>> >
>> >     OK, that got further.
>> >
>> >     On Sat, Jun 30, 2018 at 3:03 PM Mark Adams <mfadams at lbl.gov
>> >     <mailto:mfadams at lbl.gov>> wrote:
>> >
>> >         Like this?
>> >
>> >
>> >
>>  '--with-blaslapack-lib=/opt/intel/compilers_and_libraries_2018.1.163/linux/mkl/lib/intel64/libmkl_intel_thread.a',
>> >
>> >
>> >         On Sat, Jun 30, 2018 at 3:00 PM Mark Adams <mfadams at lbl.gov
>> >         <mailto:mfadams at lbl.gov>> wrote:
>> >
>> >
>> >                     Specify either "--with-blaslapack-dir" or
>> >                     "--with-blaslapack-lib --with-blaslapack-include".
>> >                     But not both!
>> >
>> >
>> >                 Get rid of the dir option, and give the full path to the
>> >                 library.
>> >
>> >
>> >             What is the syntax for giving the full path?
>> >
>>
> <configure.log>
>
>
> --
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/20180701/d33abad1/attachment.html>


More information about the petsc-dev mailing list