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,e; PetscInt m,n,j,imin[M],imax[M],imaxabs[M],indices[N],row,testcase=0; PetscScalar values[N]; PetscErrorCode ierr; PetscMPIInt size,rank; PetscReal enorm; 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 PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"max+min > PETSC_MACHINE_EPSILON "); for (j = 0; j < n; j++) { if (imin[j] != imax[j]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"imin[%D] %D != imax %D",j,imin[j],imax[j]); } /* MatGetRowMaxAbs() */ ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMaxAbs\n");CHKERRQ(ierr); ierr = MatGetRowMaxAbs(A,maxabs,NULL);CHKERRQ(ierr); ierr = MatGetRowMaxAbs(A,maxabs,imaxabs);CHKERRQ(ierr); ierr = VecView(maxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = PetscIntView(n,imaxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); /* MatGetRowMinAbs() */ ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMinAbs\n");CHKERRQ(ierr); ierr = MatGetRowMinAbs(A,min,NULL);CHKERRQ(ierr); ierr = MatGetRowMinAbs(A,min,imin);CHKERRQ(ierr); ierr = VecView(min,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = PetscIntView(n,imin,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); if (size == 1) { /* Test MatGetRowMax, MatGetRowMin and MatGetRowMaxAbs for SeqDense and MPIBAIJ matrix */ Mat Adense; Vec max_d,maxabs_d; ierr = VecDuplicate(min,&max_d);CHKERRQ(ierr); ierr = VecDuplicate(min,&maxabs_d);CHKERRQ(ierr); ierr = MatScale(A,-1.0);CHKERRQ(ierr); ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMax for seqdense matrix\n");CHKERRQ(ierr); ierr = MatGetRowMax(Adense,max_d,imax);CHKERRQ(ierr); ierr = VecWAXPY(e,-1.0,max,max_d);CHKERRQ(ierr); /* e = -max + max_d */ ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr); if (enorm > PETSC_MACHINE_EPSILON) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-max + max_d) %g > PETSC_MACHINE_EPSILON",(double)enorm); ierr = MatScale(Adense,-1.0);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMin for seqdense matrix\n");CHKERRQ(ierr); ierr = MatGetRowMin(Adense,min,imin);CHKERRQ(ierr); ierr = VecWAXPY(e,1.0,max,min);CHKERRQ(ierr); /* e = max + min */ ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr); if (enorm > PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"max+min > PETSC_MACHINE_EPSILON "); for (j = 0; j < n; j++) { if (imin[j] != imax[j]) { SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"imin[%D] %D != imax %D for seqdense matrix",j,imin[j],imax[j]); } } ierr = PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMaxAbs for seqdense matrix\n");CHKERRQ(ierr); ierr = MatGetRowMaxAbs(Adense,maxabs_d,imaxabs);CHKERRQ(ierr); ierr = VecWAXPY(e,-1.0,maxabs,maxabs_d);CHKERRQ(ierr); /* e = -maxabs + maxabs_d */ ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr); if (enorm > PETSC_MACHINE_EPSILON) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON",(double)enorm); ierr = MatDestroy(&Adense);CHKERRQ(ierr); ierr = VecDestroy(&max_d);CHKERRQ(ierr); ierr = VecDestroy(&maxabs_d);CHKERRQ(ierr); } { /* BAIJ matrix */ Mat B; Vec maxabsB,maxabsB2; PetscInt bs=2,*imaxabsB,*imaxabsB2,rstart,rend,cstart,cend,ncols,col,Brows[2],Bcols[2]; const PetscInt *cols; const PetscScalar *vals,*vals2; PetscScalar Bvals[4]; ierr = PetscMalloc2(M,&imaxabsB,bs*M,&imaxabsB2);CHKERRQ(ierr); /* bs = 1 */ ierr = MatConvert(A,MATMPIBAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); ierr = VecDuplicate(min,&maxabsB);CHKERRQ(ierr); ierr = MatGetRowMaxAbs(B,maxabsB,NULL);CHKERRQ(ierr); ierr = MatGetRowMaxAbs(B,maxabsB,imaxabsB);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMaxAbs for BAIJ matrix\n");CHKERRQ(ierr); ierr = VecWAXPY(e,-1.0,maxabs,maxabsB);CHKERRQ(ierr); /* e = -maxabs + maxabsB */ ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr); if (enorm > PETSC_MACHINE_EPSILON) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON",(double)enorm); for (j = 0; j < n; j++) { if (imaxabs[j] != imaxabsB[j]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"imaxabs[%D] %D != imaxabsB %D",j,imin[j],imax[j]); } ierr = MatDestroy(&B);CHKERRQ(ierr); /* Test bs = 2: Create B with bs*bs block structure of A */ ierr = VecCreate(PETSC_COMM_WORLD,&maxabsB2);CHKERRQ(ierr); ierr = VecSetSizes(maxabsB2,bs*m,PETSC_DECIDE);CHKERRQ(ierr); ierr = VecSetFromOptions(maxabsB2);CHKERRQ(ierr); ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); ierr = MatSetSizes(B,bs*(rend-rstart),bs*(cend-cstart),PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); ierr = MatSetFromOptions(B);CHKERRQ(ierr); ierr = MatSetUp(B);CHKERRQ(ierr); for (row=rstart; row PETSC_MACHINE_EPSILON) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row %D maxabsB != maxabsB2",row); } ierr = VecRestoreArrayRead(maxabsB,&vals);CHKERRQ(ierr); ierr = VecRestoreArrayRead(maxabsB2,&vals2);CHKERRQ(ierr); for (col=0; col