[petsc-users] MatGetLocalToGlobalMapping returning null pointers

Klinvex, Alicia Marie amklinv at sandia.gov
Mon Jul 27 14:29:50 CDT 2015


Hello,


I am trying to determine how the rows and columns of my PETSc Mat are being mapped to different MPI processes.  I call MatGetLocalToGlobalMapping, but the ISLocalToGlobalMapping obtained for both the rows and columns is null, and I receive the following error message when I call ISLocalToGlobalMappingGetSize:


[1]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
[1]PETSC ERROR: Null argument, when expecting valid pointer
[1]PETSC ERROR: Null Object: Parameter # 1
[1]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
[1]PETSC ERROR: Petsc Release Version 3.6.0, Jun, 09, 2015
[1]PETSC ERROR: ./AnasaziTest.exe on a arch-linux2-cxx-debug named hardin.srn.sandia.gov by amklinv Mon Jul 27 10:25:35 2015
[1]PETSC ERROR: Configure options --with-blas-lib=/home/amklinv/lapack-3.5.0/librefblas.a --with-lapack-lib=/home/amklinv/lapack-3.5.0/liblapack.a --with-mpi-dir=/projects/install/rhel6-x86_64/sems/compiler/gcc/4.7.2/openmpi/1.6.5 --prefix=/home/amklinv/PETSc/petsc_install --with-shared-libraries=0 --with-cxx-dialect=C++11 --with-c++-support --with-clanguage=cxx --download-superlu_dist --download-superlu --download-hypre --download-metis --download-parmetis
[1]PETSC ERROR: #1 ISLocalToGlobalMappingGetSize() line 55 in /home/amklinv/PETSc/petsc-3.6.0/src/vec/is/utils/isltog.c
[1]PETSC ERROR: #2 Tpetra_PETScAIJRowGraph() line 247 in Tpetra_PETScAIJRowGraph.hpp
[2]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------


Does this mean I set up the matrix incorrectly, or is there a different function I should be calling instead?


The relevant pieces of source code are below (although I would be happy to provide additional information if it would prove useful).


Thank you,

Alicia Klinvex




  Mat            A;        /* PETSc matrix */

  PetscInitialize(&argc,&args,(char *)0,help);

  ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
  ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
  ierr = MatSetType(A, MATAIJ);CHKERRQ(ierr);
  ierr = MatSetFromOptions(A);CHKERRQ(ierr);
  ierr = MatMPIAIJSetPreallocation(A,5,PETSC_NULL,5,PETSC_NULL);CHKERRQ(ierr);
  ierr = MatSetUp(A);CHKERRQ(ierr);
  PetscObjectGetComm( (PetscObject)A, &comm);
  ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
  if (!rank) printf("Matrix has %d (%dx%d) rows\n",m*n,m,n);

  ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);

  for (Ii=Istart; Ii<Iend; Ii++) {
    v = -1.0; i = Ii/n; j = Ii - i*n;
    if (i>0)   {J = Ii - n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
    if (i<m-1) {J = Ii + n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
    if (j>0)   {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
    if (j<n-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
    v = 4.0; ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
  }

  ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);


...



 ISLocalToGlobalMapping PETScRowOwnership, PETScColOwnership;

  ierr = MatGetLocalToGlobalMapping(A, &PETScRowOwnership, &PETScColOwnership); CHKERRV(ierr);
  ierr = ISLocalToGlobalMappingGetSize(PETScRowOwnership, &PETScLocalCols); CHKERRV(ierr);
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20150727/41404ca2/attachment.html>


More information about the petsc-users mailing list