1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatGetRowMax(), MatGetRowMin(), MatGetRowMaxAbs()\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown #define M 5 7c4762a1bSJed Brown #define N 6 8c4762a1bSJed Brown 9c4762a1bSJed Brown int main(int argc,char **args) 10c4762a1bSJed Brown { 11c4762a1bSJed Brown Mat A; 12fa213d2fSHong Zhang Vec min,max,maxabs,e; 134e879edeSHong Zhang PetscInt m,n,j,imin[M],imax[M],imaxabs[M],indices[N],row,testcase=0; 14c4762a1bSJed Brown PetscScalar values[N]; 15c4762a1bSJed Brown PetscErrorCode ierr; 161a254869SHong Zhang PetscMPIInt size,rank; 17fa213d2fSHong Zhang PetscReal enorm; 18c4762a1bSJed Brown 19c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 20ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 21ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 221a254869SHong Zhang ierr = PetscOptionsGetInt(NULL,NULL,"-testcase",&testcase,NULL);CHKERRQ(ierr); 23c4762a1bSJed Brown 24c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 251a254869SHong Zhang if (testcase == 1) { /* proc[0] holds entire A and other processes have no entry */ 26dd400576SPatrick Sanan if (rank == 0) { 271a254869SHong Zhang ierr = MatSetSizes(A,M,N,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 281a254869SHong Zhang } else { 291a254869SHong Zhang ierr = MatSetSizes(A,0,0,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 301a254869SHong Zhang } 311a254869SHong Zhang testcase = 0; 321a254869SHong Zhang } else { 33c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 341a254869SHong Zhang } 35c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 371a254869SHong Zhang 38dd400576SPatrick Sanan if (rank == 0) { /* proc[0] sets matrix A */ 391a254869SHong Zhang for (j=0; j<N; j++) indices[j] = j; 401a254869SHong Zhang switch (testcase) { 411a254869SHong Zhang case 1: /* see testcast 0 */ 421a254869SHong Zhang break; 431a254869SHong Zhang case 2: 44c4762a1bSJed Brown row = 0; 45f07e67edSHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; values[3] = -4.0; values[4] = 1.0; values[5] = 1.0; 461a254869SHong Zhang ierr = MatSetValues(A,1,&row,N,indices,values,INSERT_VALUES);CHKERRQ(ierr); 471a254869SHong Zhang row = 2; 481a254869SHong Zhang indices[0] = 0; indices[1] = 3; indices[2] = 5; 491a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 501a254869SHong Zhang ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr); 511a254869SHong Zhang row = 3; 521a254869SHong Zhang indices[0] = 0; indices[1] = 1; indices[2] = 4; 531a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 54c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr); 55c4762a1bSJed Brown row = 4; 561a254869SHong Zhang indices[0] = 0; indices[1] = 1; indices[2] = 2; 571a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 581a254869SHong Zhang ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr); 591a254869SHong Zhang break; 601a254869SHong Zhang case 3: 611a254869SHong Zhang row = 0; 621a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 631a254869SHong Zhang ierr = MatSetValues(A,1,&row,3,indices+1,values,INSERT_VALUES);CHKERRQ(ierr); 641a254869SHong Zhang row = 1; 651a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 661a254869SHong Zhang ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr); 671a254869SHong Zhang row = 2; 681a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 691a254869SHong Zhang ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr); 701a254869SHong Zhang row = 3; 711a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; values[3] = -1.0; 721a254869SHong Zhang ierr = MatSetValues(A,1,&row,4,indices,values,INSERT_VALUES);CHKERRQ(ierr); 731a254869SHong Zhang row = 4; 741a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; values[3] = -1.0; 751a254869SHong Zhang ierr = MatSetValues(A,1,&row,4,indices,values,INSERT_VALUES);CHKERRQ(ierr); 761a254869SHong Zhang break; 771a254869SHong Zhang 781a254869SHong Zhang default: 791a254869SHong Zhang row = 0; 801a254869SHong Zhang values[0] = -1.0; values[1] = 0.0; values[2] = 1.0; values[3] = 3.0; values[4] = 4.0; values[5] = -5.0; 811a254869SHong Zhang ierr = MatSetValues(A,1,&row,N,indices,values,INSERT_VALUES);CHKERRQ(ierr); 821a254869SHong Zhang row = 1; 831a254869SHong Zhang ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr); 841a254869SHong Zhang row = 3; 85c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,indices+4,values+4,INSERT_VALUES);CHKERRQ(ierr); 86c4762a1bSJed Brown row = 4; 87c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,2,indices+4,values+4,INSERT_VALUES);CHKERRQ(ierr); 881a254869SHong Zhang } 891a254869SHong Zhang } 90c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 93c4762a1bSJed Brown 94c4762a1bSJed Brown ierr = MatGetLocalSize(A, &m,&n);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&min);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = VecSetSizes(min,m,PETSC_DECIDE);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = VecSetFromOptions(min);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = VecDuplicate(min,&max);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = VecDuplicate(min,&maxabs);CHKERRQ(ierr); 100fa213d2fSHong Zhang ierr = VecDuplicate(min,&e);CHKERRQ(ierr); 101c4762a1bSJed Brown 102f07e67edSHong Zhang /* MatGetRowMax() */ 103f07e67edSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMax\n");CHKERRQ(ierr); 1041a254869SHong Zhang ierr = MatGetRowMax(A,max,NULL);CHKERRQ(ierr); 1051a254869SHong Zhang ierr = MatGetRowMax(A,max,imax);CHKERRQ(ierr); 1061a254869SHong Zhang ierr = VecView(max,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1071a254869SHong Zhang ierr = VecGetLocalSize(max,&n);CHKERRQ(ierr); 1081a254869SHong Zhang ierr = PetscIntView(n,imax,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 109fa213d2fSHong Zhang 110f07e67edSHong Zhang /* MatGetRowMin() */ 111fa213d2fSHong Zhang ierr = MatScale(A,-1.0);CHKERRQ(ierr); 112f07e67edSHong Zhang ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 113f07e67edSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMin\n");CHKERRQ(ierr); 114fa213d2fSHong Zhang ierr = MatGetRowMin(A,min,NULL);CHKERRQ(ierr); 115fa213d2fSHong Zhang ierr = MatGetRowMin(A,min,imin);CHKERRQ(ierr); 116fa213d2fSHong Zhang 117fa213d2fSHong Zhang ierr = VecWAXPY(e,1.0,max,min);CHKERRQ(ierr); /* e = max + min */ 118fa213d2fSHong Zhang ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr); 119*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm > PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"max+min > PETSC_MACHINE_EPSILON "); 120f07e67edSHong Zhang for (j = 0; j < n; j++) { 121*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(imin[j] != imax[j],PETSC_COMM_SELF,PETSC_ERR_PLIB,"imin[%" PetscInt_FMT "] %" PetscInt_FMT " != imax %" PetscInt_FMT,j,imin[j],imax[j]); 122f07e67edSHong Zhang } 123fa213d2fSHong Zhang 124f07e67edSHong Zhang /* MatGetRowMaxAbs() */ 125f07e67edSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMaxAbs\n");CHKERRQ(ierr); 126475b8b61SHong Zhang ierr = MatGetRowMaxAbs(A,maxabs,NULL);CHKERRQ(ierr); 127475b8b61SHong Zhang ierr = MatGetRowMaxAbs(A,maxabs,imaxabs);CHKERRQ(ierr); 128475b8b61SHong Zhang ierr = VecView(maxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 129475b8b61SHong Zhang ierr = PetscIntView(n,imaxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1301a254869SHong Zhang 131f07e67edSHong Zhang /* MatGetRowMinAbs() */ 132f07e67edSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMinAbs\n");CHKERRQ(ierr); 133f07e67edSHong Zhang ierr = MatGetRowMinAbs(A,min,NULL);CHKERRQ(ierr); 134f07e67edSHong Zhang ierr = MatGetRowMinAbs(A,min,imin);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = VecView(min,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 136f07e67edSHong Zhang ierr = PetscIntView(n,imin,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 137f07e67edSHong Zhang 138f07e67edSHong Zhang if (size == 1) { 139f07e67edSHong Zhang /* Test MatGetRowMax, MatGetRowMin and MatGetRowMaxAbs for SeqDense and MPIBAIJ matrix */ 1404e879edeSHong Zhang Mat Adense; 1414e879edeSHong Zhang Vec max_d,maxabs_d; 1424e879edeSHong Zhang ierr = VecDuplicate(min,&max_d);CHKERRQ(ierr); 1434e879edeSHong Zhang ierr = VecDuplicate(min,&maxabs_d);CHKERRQ(ierr); 144f07e67edSHong Zhang 145f07e67edSHong Zhang ierr = MatScale(A,-1.0);CHKERRQ(ierr); 1464e879edeSHong Zhang ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);CHKERRQ(ierr); 1474e879edeSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMax for seqdense matrix\n");CHKERRQ(ierr); 1484e879edeSHong Zhang ierr = MatGetRowMax(Adense,max_d,imax);CHKERRQ(ierr); 1494e879edeSHong Zhang 1504e879edeSHong Zhang ierr = VecWAXPY(e,-1.0,max,max_d);CHKERRQ(ierr); /* e = -max + max_d */ 1514e879edeSHong Zhang ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr); 152*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm > PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-max + max_d) %g > PETSC_MACHINE_EPSILON",(double)enorm); 1534e879edeSHong Zhang 1544e879edeSHong Zhang ierr = MatScale(Adense,-1.0);CHKERRQ(ierr); 1554e879edeSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMin for seqdense matrix\n");CHKERRQ(ierr); 1564e879edeSHong Zhang ierr = MatGetRowMin(Adense,min,imin);CHKERRQ(ierr); 157f07e67edSHong Zhang 158f07e67edSHong Zhang ierr = VecWAXPY(e,1.0,max,min);CHKERRQ(ierr); /* e = max + min */ 159f07e67edSHong Zhang ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr); 160*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm > PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"max+min > PETSC_MACHINE_EPSILON "); 161f07e67edSHong Zhang for (j = 0; j < n; j++) { 162f07e67edSHong Zhang if (imin[j] != imax[j]) { 16398921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"imin[%" PetscInt_FMT "] %" PetscInt_FMT " != imax %" PetscInt_FMT " for seqdense matrix",j,imin[j],imax[j]); 164f07e67edSHong Zhang } 165f07e67edSHong Zhang } 166f07e67edSHong Zhang 1674e879edeSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMaxAbs for seqdense matrix\n");CHKERRQ(ierr); 1684e879edeSHong Zhang ierr = MatGetRowMaxAbs(Adense,maxabs_d,imaxabs);CHKERRQ(ierr); 1694e879edeSHong Zhang ierr = VecWAXPY(e,-1.0,maxabs,maxabs_d);CHKERRQ(ierr); /* e = -maxabs + maxabs_d */ 1704e879edeSHong Zhang ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr); 171*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm > PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON",(double)enorm); 172c4762a1bSJed Brown 1734e879edeSHong Zhang ierr = MatDestroy(&Adense);CHKERRQ(ierr); 1744e879edeSHong Zhang ierr = VecDestroy(&max_d);CHKERRQ(ierr); 1754e879edeSHong Zhang ierr = VecDestroy(&maxabs_d);CHKERRQ(ierr); 17643359b5eSHong Zhang } 17743359b5eSHong Zhang 17843359b5eSHong Zhang { /* BAIJ matrix */ 1794e879edeSHong Zhang Mat B; 18043359b5eSHong Zhang Vec maxabsB,maxabsB2; 18143359b5eSHong Zhang PetscInt bs=2,*imaxabsB,*imaxabsB2,rstart,rend,cstart,cend,ncols,col,Brows[2],Bcols[2]; 18243359b5eSHong Zhang const PetscInt *cols; 18343359b5eSHong Zhang const PetscScalar *vals,*vals2; 18443359b5eSHong Zhang PetscScalar Bvals[4]; 18543359b5eSHong Zhang 18643359b5eSHong Zhang ierr = PetscMalloc2(M,&imaxabsB,bs*M,&imaxabsB2);CHKERRQ(ierr); 18743359b5eSHong Zhang 18843359b5eSHong Zhang /* bs = 1 */ 1894e879edeSHong Zhang ierr = MatConvert(A,MATMPIBAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 1904e879edeSHong Zhang ierr = VecDuplicate(min,&maxabsB);CHKERRQ(ierr); 19143359b5eSHong Zhang ierr = MatGetRowMaxAbs(B,maxabsB,NULL);CHKERRQ(ierr); 19243359b5eSHong Zhang ierr = MatGetRowMaxAbs(B,maxabsB,imaxabsB);CHKERRQ(ierr); 19343359b5eSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMaxAbs for BAIJ matrix\n");CHKERRQ(ierr); 1944e879edeSHong Zhang ierr = VecWAXPY(e,-1.0,maxabs,maxabsB);CHKERRQ(ierr); /* e = -maxabs + maxabsB */ 1954e879edeSHong Zhang ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr); 196*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm > PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON",(double)enorm); 1974e879edeSHong Zhang 1984e879edeSHong Zhang for (j = 0; j < n; j++) { 199*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(imaxabs[j] != imaxabsB[j],PETSC_COMM_SELF,PETSC_ERR_PLIB,"imaxabs[%" PetscInt_FMT "] %" PetscInt_FMT " != imaxabsB %" PetscInt_FMT,j,imin[j],imax[j]); 200c4762a1bSJed Brown } 2014e879edeSHong Zhang ierr = MatDestroy(&B);CHKERRQ(ierr); 20243359b5eSHong Zhang 20343359b5eSHong Zhang /* Test bs = 2: Create B with bs*bs block structure of A */ 20443359b5eSHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&maxabsB2);CHKERRQ(ierr); 20543359b5eSHong Zhang ierr = VecSetSizes(maxabsB2,bs*m,PETSC_DECIDE);CHKERRQ(ierr); 20643359b5eSHong Zhang ierr = VecSetFromOptions(maxabsB2);CHKERRQ(ierr); 20743359b5eSHong Zhang 20843359b5eSHong Zhang ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 20943359b5eSHong Zhang ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr); 21043359b5eSHong Zhang ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 21143359b5eSHong Zhang ierr = MatSetSizes(B,bs*(rend-rstart),bs*(cend-cstart),PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 21243359b5eSHong Zhang ierr = MatSetFromOptions(B);CHKERRQ(ierr); 21343359b5eSHong Zhang ierr = MatSetUp(B);CHKERRQ(ierr); 21443359b5eSHong Zhang 21543359b5eSHong Zhang for (row=rstart; row<rend; row++) { 21643359b5eSHong Zhang ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 21743359b5eSHong Zhang for (col=0; col<ncols; col++) { 21843359b5eSHong Zhang for (j=0; j<bs; j++) { 21943359b5eSHong Zhang Brows[j] = bs*row + j; 22043359b5eSHong Zhang Bcols[j] = bs*cols[col]+j; 22143359b5eSHong Zhang } 22243359b5eSHong Zhang for (j=0; j<bs*bs; j++) Bvals[j] = vals[col]; 22343359b5eSHong Zhang ierr = MatSetValues(B,bs,Brows,bs,Bcols,Bvals,INSERT_VALUES);CHKERRQ(ierr); 22443359b5eSHong Zhang } 22543359b5eSHong Zhang ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 22643359b5eSHong Zhang } 22743359b5eSHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 22843359b5eSHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 22943359b5eSHong Zhang 23043359b5eSHong Zhang ierr = MatGetRowMaxAbs(B,maxabsB2,imaxabsB2);CHKERRQ(ierr); 23143359b5eSHong Zhang 23243359b5eSHong Zhang /* Check maxabsB2 and imaxabsB2 */ 23343359b5eSHong Zhang ierr = VecGetArrayRead(maxabsB,&vals);CHKERRQ(ierr); 23443359b5eSHong Zhang ierr = VecGetArrayRead(maxabsB2,&vals2);CHKERRQ(ierr); 23543359b5eSHong Zhang for (row=0; row<m; row++) { 23643359b5eSHong Zhang if (PetscAbsScalar(vals[row] - vals2[bs*row]) > PETSC_MACHINE_EPSILON) 23798921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row %" PetscInt_FMT " maxabsB != maxabsB2",row); 23843359b5eSHong Zhang } 23943359b5eSHong Zhang ierr = VecRestoreArrayRead(maxabsB,&vals);CHKERRQ(ierr); 24043359b5eSHong Zhang ierr = VecRestoreArrayRead(maxabsB2,&vals2);CHKERRQ(ierr); 24143359b5eSHong Zhang 24243359b5eSHong Zhang for (col=0; col<n; col++) { 24343359b5eSHong Zhang if (imaxabsB[col] != imaxabsB2[bs*col]/bs) 24498921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"col %" PetscInt_FMT " imaxabsB != imaxabsB2",col); 24543359b5eSHong Zhang } 2464e879edeSHong Zhang ierr = VecDestroy(&maxabsB);CHKERRQ(ierr); 24743359b5eSHong Zhang ierr = MatDestroy(&B);CHKERRQ(ierr); 24843359b5eSHong Zhang ierr = VecDestroy(&maxabsB2);CHKERRQ(ierr); 24943359b5eSHong Zhang ierr = PetscFree2(imaxabsB,imaxabsB2);CHKERRQ(ierr); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252c4762a1bSJed Brown ierr = VecDestroy(&min);CHKERRQ(ierr); 253c4762a1bSJed Brown ierr = VecDestroy(&max);CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = VecDestroy(&maxabs);CHKERRQ(ierr); 255fa213d2fSHong Zhang ierr = VecDestroy(&e);CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = PetscFinalize(); 258c4762a1bSJed Brown return ierr; 259c4762a1bSJed Brown } 260c4762a1bSJed Brown 261c4762a1bSJed Brown /*TEST 262c4762a1bSJed Brown 263c4762a1bSJed Brown test: 264c4762a1bSJed Brown output_file: output/ex114.out 265c4762a1bSJed Brown 266c4762a1bSJed Brown test: 267c4762a1bSJed Brown suffix: 2 2684e879edeSHong Zhang args: -testcase 1 2694e879edeSHong Zhang output_file: output/ex114.out 270c4762a1bSJed Brown 271c4762a1bSJed Brown test: 272c4762a1bSJed Brown suffix: 3 2734e879edeSHong Zhang args: -testcase 2 2744e879edeSHong Zhang output_file: output/ex114_3.out 2754e879edeSHong Zhang 2764e879edeSHong Zhang test: 2774e879edeSHong Zhang suffix: 4 2784e879edeSHong Zhang args: -testcase 3 2794e879edeSHong Zhang output_file: output/ex114_4.out 2804e879edeSHong Zhang 2814e879edeSHong Zhang test: 2824e879edeSHong Zhang suffix: 5 2834e879edeSHong Zhang nsize: 3 2844e879edeSHong Zhang args: -testcase 0 2854e879edeSHong Zhang output_file: output/ex114_5.out 2864e879edeSHong Zhang 2874e879edeSHong Zhang test: 2884e879edeSHong Zhang suffix: 6 2894e879edeSHong Zhang nsize: 3 2904e879edeSHong Zhang args: -testcase 1 2914e879edeSHong Zhang output_file: output/ex114_6.out 2924e879edeSHong Zhang 2934e879edeSHong Zhang test: 2944e879edeSHong Zhang suffix: 7 2954e879edeSHong Zhang nsize: 3 2964e879edeSHong Zhang args: -testcase 2 2974e879edeSHong Zhang output_file: output/ex114_7.out 2984e879edeSHong Zhang 2994e879edeSHong Zhang test: 3004e879edeSHong Zhang suffix: 8 3014e879edeSHong Zhang nsize: 3 3024e879edeSHong Zhang args: -testcase 3 3034e879edeSHong Zhang output_file: output/ex114_8.out 304c4762a1bSJed Brown 305c4762a1bSJed Brown TEST*/ 306