[petsc-users] Loading a PETsc matrix with matrix data in CSR format?

Barry Smith bsmith at mcs.anl.gov
Tue Oct 3 08:07:09 CDT 2017


> On Oct 3, 2017, at 1:40 PM, Klaus Burkart <k_burkart at yahoo.com> wrote:
> 
> Hello, 
> 
> I am still struggling to load my matrix in a PETSc matrix.
> 
> I tried the following code which is part of a c++ function which converts the original matrix to CSR fromat and should load a PETSc matrix.
> 
>     Mat M;
>    
  Cannot have the following line

> MatSetFromOptions(M);
>     // fill PETSc matrix

Should be &M at the end not M

>     MatCreateSeqAIJWithArrays(PETSC_COMM_WORLD, n, n, rows, cols, vals, M)
> 
> but I get the following error message when compiling the code:
> 
>  error: cannot convert ‘Mat {aka _p_Mat*}’ to ‘_p_Mat**’ for argument ‘7’ to ‘PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*, _p_Mat**)’
>      MatCreateSeqAIJWithArrays(PETSC_COMM_WORLD, n, n, rows, cols, vals, M)
> 
> What's wrong with this?
> 
> Klaus
> 
> 
> Barry Smith <bsmith at mcs.anl.gov> schrieb am 16:38 Sonntag, 17.September 2017:
> 
> 
> 
> 1)  MatSetValues() has absolutely nothing to do with the format used internally by PETSc to store the matrix. MatSetValues() is used to convey blocks of values into the matrix. 
> 
>   2) If you want to load a matrix onto multiple processes from a file NEVER write your own parallel ASCII matrix reader. Instead in either C, C++, Python, or Matlab write a routine that reads in the matrix SEQUENTIALLY and then saves it with MatView() and a binary viewer. You can then load the matrix easily and efficiently in PETSc in parallel with MatLoad
> 
>   3) If you have matrix already in CSR format you can use MatCreateSeqAIJWithArrays()
> 
>   Barry
> 
> 
> 
> > On Sep 17, 2017, at 9:25 AM, Klaus Burkart <k_burkart at yahoo.com> wrote:
> > 
> > The matrix import function looks like this:
> > 
> > void csr2pet
> > (
> >    const Foam::lduMatrix & matrix,
> >    petsc_declaration<Type> & petsc_matrix  // How to declare the PETsc matrix to be filled?
> > )
> > {
> >    int n = matrix.diag().size(); // small case n = 40800
> >    int nnz = matrix.lower().size() + matrix.upper().size() + matrix.diag().size(); // small case nnz = 203800
> > 
> >    // allocate memory for CSR sparse matrix using calloc
> >    ScalarType * vals = (ScalarType *)calloc(nnz, sizeof(ScalarType));
> >    uint * cols = (uint *)calloc(nnz, sizeof(uint));
> >    uint * rows = (uint *)calloc(n, sizeof(uint));
> > 
> >    // call function to convert original LDU matrix to CSR format
> >    exPet::ldu2csr(matrix,rows,cols,vals);
> > 
> >    // fill PETsc matrix
> >    MatSetValues(petsc_matrix, ?, ?, ?, ?, ?, INSERT_VALUES);
> > 
> >    // free and release the matrix memory
> >    free(rows); free(cols); free(vals);  // calloc()
> > }
> > 
> > 
> > Questions:
> > 
> > 1: How to declare the petsc_matrix to be filled by the function with the content of the original matrix?
> > 
> > 2: MatSetValues(petsc_matrix, ?, ?, ?, ?, ?, INSERT_VALUES); is used to actually fill the petsc_matrix and I was of the opinion that PETsc uses the CSR format but I can't work out how CSR format is described by:
> > 
> >    v        - a logically two-dimensional array of values
> >    m, idxm    - the number of rows and their global indices
> >    n, idxn    - the number of columns and their global indices
> > 
> > My original matrix is converted to CSR format, i.e. three arrays cols (column_indices), rows (row_start_indices) and vals (values).
> > 
> > How can I load my matrix into a PETsc matrix for parallel processing? MatSetValues(petsc_matrix, ?, ?, ?, ?, ?, INSERT_VALUES);
> > 
> > Klaus
> 
> 



More information about the petsc-users mailing list