[petsc-users] mkl_cpardiso: Matrix B must be MATSEQDENSE matrix
Smith, Barry F.
bsmith at mcs.anl.gov
Fri Aug 31 02:36:01 CDT 2018
Marius,
Try editing src/mat/impls/aij/mpi/mkl_cpardiso/mkl_cpardiso.c at the subroutine MatMatSolve_MKL_CPARDISO() and remove the lines
ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&flg);CHKERRQ(ierr);
if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATSEQDENSE matrix");
ierr = PetscObjectTypeCompare((PetscObject)X,MATSEQDENSE,&flg);CHKERRQ(ierr);
if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATSEQDENSE matrix");
Run
make gnumake
then try compiling and running your application again.
Barry
It could be that this code was cut and pasted from some other solver (Mumps) and hence has checks for SEQDENSE that are not needed for MKL_CPARDISO.
> On Aug 31, 2018, at 2:03 AM, Marius Buerkle <mbuerkle at web.de> wrote:
>
> Hi !
>
> I tried to use MKL_CPARDISO solver in parallel but I get the following error
>
> [0]PETSC ERROR: Invalid argument
> [0]PETSC ERROR: Matrix B must be MATSEQDENSE matrix
> [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.3-1202-ge6d955127a GIT Date: 2018-08-23 16:59:59 -0500
> [0]PETSC ERROR: /home/marius/prog/ownstuff/fortran/programs/transomat_dev/transomat/transomat on a named tono-hpc1 by marius Fri Aug 31 15:57:15 2018
> [0]PETSC ERROR: Configure options --prefix=/home/marius/prog/petsc/master_opt --with-64-bit-indices=0 --CC=mpicc --COPTFLAGS="-O3 -march=native -g -std=c11" --CXX=mpicxx --CXXOPTFLAGS="-O3 -march=native -g -std=c++11" --FC=mpiifort --FOPTFLAGS="-O2 -xHost -g" --with-mpi=1 --with-x=1 --download-parmetis=yes --download-metis=yes --download-ptscotch=yes --download-scotch=yes --download-mumps=yes --with-blaslapack-lib="-L/home/marius/intel/compilers_and_libraries_2018.3.222/linux/mkl/lib/intel64 -lmkl_scalapack_lp64 -lmkl_intel_lp64 -lmkl_intel_thread -lmkl_core -lmkl_blacs_intelmpi_lp64 -liomp5 -lpthread -lm -ldl" --with-scalapack-lib="-L/home/marius/intel/compilers_and_libraries_2018.3.222/linux/mkl/lib/intel64 -lmkl_scalapack_lp64 -lmkl_intel_lp64 -lmkl_intel_thread -lmkl_core -lmkl_blacs_intelmpi_lp64 -liomp5 -lpthread -lm -ldl" --download-pastix=no --download-superlu_dist=yes --download-superlu_dist-commit=803c8ab665b33e9136b199c4521be465e026bb0d --download-hwloc=yes --with-mkl_pardiso-dir=/home/marius/intel/compilers_and_libraries_2018.3.222/linux/mkl --with-openmp=0 --with-pthread=0 --download-elemental=yes --download-elemental-commit=de7b5bea1abf5f626b91582f742cf99e2e551bff --with-mkl_cpardiso-dir=/home/marius/intel/compilers_and_libraries_2018.3.222/linux/mkl --with-mkl_sparse_optimize-dir=/home/marius/intel/compilers_and_libraries_2018.3.222/linux/mkl --with-mkl_sparse_sp2m-dir=/home/marius/intel/compilers_and_libraries_2018.3.222/linux/mkl --download-sowing=yes --with-cxx-dialect=c++11 --with-scalar-type=complex --with-debugging=0
> [0]PETSC ERROR: #1 MatMatSolve_MKL_CPARDISO() line 319 in /home/marius/prog/petsc/petsc/src/mat/impls/aij/mpi/mkl_cpardiso/mkl_cpardiso.c
> [0]PETSC ERROR: #2 MatMatSolve() line 3385 in /home/marius/prog/petsc/petsc/src/mat/interface/matrix.c
>
>
> I wonder isn't CPARDISO supposed to work with MATMPIDENSE ?
>
> best,
> marius
More information about the petsc-users
mailing list