static char help[] = "Tests MatGetRowMax(), MatGetRowMin(), MatGetRowMaxAbs()\n"; #include #define M 5 #define N 6 int main(int argc,char **args) { Mat A; Vec min,max,maxabs; PetscInt m,n; PetscInt imin[M],imax[M],imaxabs[M],indices[N],row,testcase=0; PetscScalar values[N]; PetscErrorCode ierr; MatType type; PetscMPIInt size,rank; PetscBool doTest=PETSC_TRUE; PetscInt j,col,cols[N]; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-testcase",&testcase,NULL);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); if (testcase == 1) { /* proc[0] holds entire A and other processes have no entry */ if (!rank) { ierr = MatSetSizes(A,M,N,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); } else { ierr = MatSetSizes(A,0,0,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); } testcase = 0; } else { ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); } ierr = MatSetFromOptions(A);CHKERRQ(ierr); ierr = MatSetUp(A);CHKERRQ(ierr); if (!rank) { /* proc[0] sets matrix A */ for (j=0; j