[petsc-users] Extracting data from a Petsc matrix

Barry Smith bsmith at mcs.anl.gov
Sat Jul 13 11:54:00 CDT 2013


   If MAGMA "implements LAPACK functions on a CPU+GPU based system" using the standard BLAS/LAPACK interface then you should be able to simply ./configure PETSc to use the MAGMA BLAS/LAPACK library instead of the standard BLAS/LAPACK library see http://www.mcs.anl.gov/petsc/documentation/installation.html#blas-lapack  because PETSc for SeqDense matrices simply calls the BLAS/LAPACK routines, see for example MatLUFactor_SeqDense() which simply calls 

  ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
  ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
  if (!mat->pivots) {
    ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
    ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
  }
  if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
  ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
  PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info));
  ierr = PetscFPTrapPop();CHKERRQ(ierr);

  if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
  if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");

If MAGMA actually changes the calling sequences of the BLAS/LAPACK routines then it won't work automatically, though it may be possible to modify the routines in PETSc's dense.c source code to call the MAGMA versions instead of the "standard" BLAS/LAPACK routines.

    Barry

On Jul 13, 2013, at 11:43 AM, Jed Brown <jedbrown at mcs.anl.gov> wrote:

> "hsahasra at purdue.edu" <hsahasra at purdue.edu> writes:
> 
>> Hi,
>> 
>> I am working on solving a system of linear equations with square
>> matrix. I'm first factoring the matrix using LU decomposition.
> 
> I assume you're solving a dense problem because that is all MAGMA does.
> 
>> I want to do the LU decomposition step using MAGMA on GPUs. MAGMA
>> library implements LAPACK functions on a CPU+GPU based system.
>> 
>> So my question is, how do I extract the data from a Petsc Mat so that
>> it can be sent to the dgetrf routine in MAGMA.
> 
> MatDenseGetArray
> 
>> Is there any need for duplicating the data for this step?
> 
> You're on your own for storage of factors.  Alternatively, you could add
> library support so that you could use PCLU and
> '-pc_factor_mat_solver_package magma' (or PCFactorSetMatSolverPackage).
> Doing this is not a priority for us, but we can provide guidance if you
> want to tackle it.



More information about the petsc-users mailing list