static char help[] = "Test conversion of ScaLAPACK matrices.\n\n"; #include int main(int argc, char** argv) { Mat A,A_scalapack; PetscInt i,j,M=10,N=5,nloc,mloc,nrows,ncols; PetscErrorCode ierr; PetscMPIInt rank,size; IS isrows,iscols; const PetscInt *rows,*cols; PetscScalar *v; MatType type; PetscBool isDense,isAIJ,flg; ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); /* Create a matrix */ ierr = MatCreate(PETSC_COMM_WORLD, &A);CHKERRQ(ierr); mloc = PETSC_DECIDE; ierr = PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&mloc,&M);CHKERRQ(ierr); nloc = PETSC_DECIDE; ierr = PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&nloc,&N);CHKERRQ(ierr); ierr = MatSetSizes(A,mloc,nloc,M,N);CHKERRQ(ierr); ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); ierr = MatSetFromOptions(A);CHKERRQ(ierr); ierr = MatSetUp(A);CHKERRQ(ierr); /* Set local matrix entries */ ierr = MatGetOwnershipIS(A,&isrows,&iscols);CHKERRQ(ierr); ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr); ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr); ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr); ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr); for (i=0; i