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]; 151a254869SHong Zhang PetscMPIInt size,rank; 16fa213d2fSHong Zhang PetscReal enorm; 17c4762a1bSJed Brown 18*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 195f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-testcase",&testcase,NULL)); 22c4762a1bSJed Brown 235f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 241a254869SHong Zhang if (testcase == 1) { /* proc[0] holds entire A and other processes have no entry */ 25dd400576SPatrick Sanan if (rank == 0) { 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,M,N,PETSC_DECIDE,PETSC_DECIDE)); 271a254869SHong Zhang } else { 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,0,0,PETSC_DECIDE,PETSC_DECIDE)); 291a254869SHong Zhang } 301a254869SHong Zhang testcase = 0; 311a254869SHong Zhang } else { 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 331a254869SHong Zhang } 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 361a254869SHong Zhang 37dd400576SPatrick Sanan if (rank == 0) { /* proc[0] sets matrix A */ 381a254869SHong Zhang for (j=0; j<N; j++) indices[j] = j; 391a254869SHong Zhang switch (testcase) { 401a254869SHong Zhang case 1: /* see testcast 0 */ 411a254869SHong Zhang break; 421a254869SHong Zhang case 2: 43c4762a1bSJed Brown row = 0; 44f07e67edSHong 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; 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,N,indices,values,INSERT_VALUES)); 461a254869SHong Zhang row = 2; 471a254869SHong Zhang indices[0] = 0; indices[1] = 3; indices[2] = 5; 481a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES)); 501a254869SHong Zhang row = 3; 511a254869SHong Zhang indices[0] = 0; indices[1] = 1; indices[2] = 4; 521a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES)); 54c4762a1bSJed Brown row = 4; 551a254869SHong Zhang indices[0] = 0; indices[1] = 1; indices[2] = 2; 561a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES)); 581a254869SHong Zhang break; 591a254869SHong Zhang case 3: 601a254869SHong Zhang row = 0; 611a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,3,indices+1,values,INSERT_VALUES)); 631a254869SHong Zhang row = 1; 641a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES)); 661a254869SHong Zhang row = 2; 671a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES)); 691a254869SHong Zhang row = 3; 701a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; values[3] = -1.0; 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,4,indices,values,INSERT_VALUES)); 721a254869SHong Zhang row = 4; 731a254869SHong Zhang values[0] = -2.0; values[1] = -2.0; values[2] = -2.0; values[3] = -1.0; 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,4,indices,values,INSERT_VALUES)); 751a254869SHong Zhang break; 761a254869SHong Zhang 771a254869SHong Zhang default: 781a254869SHong Zhang row = 0; 791a254869SHong 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; 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,N,indices,values,INSERT_VALUES)); 811a254869SHong Zhang row = 1; 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES)); 831a254869SHong Zhang row = 3; 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,1,indices+4,values+4,INSERT_VALUES)); 85c4762a1bSJed Brown row = 4; 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,2,indices+4,values+4,INSERT_VALUES)); 871a254869SHong Zhang } 881a254869SHong Zhang } 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 92c4762a1bSJed Brown 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A, &m,&n)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&min)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(min,m,PETSC_DECIDE)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(min)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(min,&max)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(min,&maxabs)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(min,&e)); 100c4762a1bSJed Brown 101f07e67edSHong Zhang /* MatGetRowMax() */ 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMax\n")); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMax(A,max,NULL)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMax(A,max,imax)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(max,PETSC_VIEWER_STDOUT_WORLD)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(max,&n)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscIntView(n,imax,PETSC_VIEWER_STDOUT_WORLD)); 108fa213d2fSHong Zhang 109f07e67edSHong Zhang /* MatGetRowMin() */ 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,-1.0)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMin\n")); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMin(A,min,NULL)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMin(A,min,imin)); 115fa213d2fSHong Zhang 1165f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(e,1.0,max,min)); /* e = max + min */ 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(e,NORM_INFINITY,&enorm)); 1182c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm > PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"max+min > PETSC_MACHINE_EPSILON "); 119f07e67edSHong Zhang for (j = 0; j < n; j++) { 1202c71b3e2SJacob 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]); 121f07e67edSHong Zhang } 122fa213d2fSHong Zhang 123f07e67edSHong Zhang /* MatGetRowMaxAbs() */ 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMaxAbs\n")); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(A,maxabs,NULL)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(A,maxabs,imaxabs)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(maxabs,PETSC_VIEWER_STDOUT_WORLD)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscIntView(n,imaxabs,PETSC_VIEWER_STDOUT_WORLD)); 1291a254869SHong Zhang 130f07e67edSHong Zhang /* MatGetRowMinAbs() */ 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMinAbs\n")); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMinAbs(A,min,NULL)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMinAbs(A,min,imin)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(min,PETSC_VIEWER_STDOUT_WORLD)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscIntView(n,imin,PETSC_VIEWER_STDOUT_WORLD)); 136f07e67edSHong Zhang 137f07e67edSHong Zhang if (size == 1) { 138f07e67edSHong Zhang /* Test MatGetRowMax, MatGetRowMin and MatGetRowMaxAbs for SeqDense and MPIBAIJ matrix */ 1394e879edeSHong Zhang Mat Adense; 1404e879edeSHong Zhang Vec max_d,maxabs_d; 1415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(min,&max_d)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(min,&maxabs_d)); 143f07e67edSHong Zhang 1445f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,-1.0)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMax for seqdense matrix\n")); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMax(Adense,max_d,imax)); 1484e879edeSHong Zhang 1495f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(e,-1.0,max,max_d)); /* e = -max + max_d */ 1505f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(e,NORM_INFINITY,&enorm)); 1512c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm > PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-max + max_d) %g > PETSC_MACHINE_EPSILON",(double)enorm); 1524e879edeSHong Zhang 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(Adense,-1.0)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMin for seqdense matrix\n")); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMin(Adense,min,imin)); 156f07e67edSHong Zhang 1575f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(e,1.0,max,min)); /* e = max + min */ 1585f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(e,NORM_INFINITY,&enorm)); 1592c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm > PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"max+min > PETSC_MACHINE_EPSILON "); 160f07e67edSHong Zhang for (j = 0; j < n; j++) { 161f07e67edSHong Zhang if (imin[j] != imax[j]) { 16298921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"imin[%" PetscInt_FMT "] %" PetscInt_FMT " != imax %" PetscInt_FMT " for seqdense matrix",j,imin[j],imax[j]); 163f07e67edSHong Zhang } 164f07e67edSHong Zhang } 165f07e67edSHong Zhang 1665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMaxAbs for seqdense matrix\n")); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(Adense,maxabs_d,imaxabs)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(e,-1.0,maxabs,maxabs_d)); /* e = -maxabs + maxabs_d */ 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(e,NORM_INFINITY,&enorm)); 1702c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm > PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON",(double)enorm); 171c4762a1bSJed Brown 1725f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Adense)); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&max_d)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&maxabs_d)); 17543359b5eSHong Zhang } 17643359b5eSHong Zhang 17743359b5eSHong Zhang { /* BAIJ matrix */ 1784e879edeSHong Zhang Mat B; 17943359b5eSHong Zhang Vec maxabsB,maxabsB2; 18043359b5eSHong Zhang PetscInt bs=2,*imaxabsB,*imaxabsB2,rstart,rend,cstart,cend,ncols,col,Brows[2],Bcols[2]; 18143359b5eSHong Zhang const PetscInt *cols; 18243359b5eSHong Zhang const PetscScalar *vals,*vals2; 18343359b5eSHong Zhang PetscScalar Bvals[4]; 18443359b5eSHong Zhang 1855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(M,&imaxabsB,bs*M,&imaxabsB2)); 18643359b5eSHong Zhang 18743359b5eSHong Zhang /* bs = 1 */ 1885f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATMPIBAIJ,MAT_INITIAL_MATRIX,&B)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(min,&maxabsB)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(B,maxabsB,NULL)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(B,maxabsB,imaxabsB)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMaxAbs for BAIJ matrix\n")); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(e,-1.0,maxabs,maxabsB)); /* e = -maxabs + maxabsB */ 1945f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(e,NORM_INFINITY,&enorm)); 1952c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm > PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON",(double)enorm); 1964e879edeSHong Zhang 1974e879edeSHong Zhang for (j = 0; j < n; j++) { 1982c71b3e2SJacob 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]); 199c4762a1bSJed Brown } 2005f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 20143359b5eSHong Zhang 20243359b5eSHong Zhang /* Test bs = 2: Create B with bs*bs block structure of A */ 2035f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&maxabsB2)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(maxabsB2,bs*m,PETSC_DECIDE)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(maxabsB2)); 20643359b5eSHong Zhang 2075f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRangeColumn(A,&cstart,&cend)); 2095f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,bs*(rend-rstart),bs*(cend-cstart),PETSC_DECIDE,PETSC_DECIDE)); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 21343359b5eSHong Zhang 21443359b5eSHong Zhang for (row=rstart; row<rend; row++) { 2155f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(A,row,&ncols,&cols,&vals)); 21643359b5eSHong Zhang for (col=0; col<ncols; col++) { 21743359b5eSHong Zhang for (j=0; j<bs; j++) { 21843359b5eSHong Zhang Brows[j] = bs*row + j; 21943359b5eSHong Zhang Bcols[j] = bs*cols[col]+j; 22043359b5eSHong Zhang } 22143359b5eSHong Zhang for (j=0; j<bs*bs; j++) Bvals[j] = vals[col]; 2225f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,bs,Brows,bs,Bcols,Bvals,INSERT_VALUES)); 22343359b5eSHong Zhang } 2245f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(A,row,&ncols,&cols,&vals)); 22543359b5eSHong Zhang } 2265f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 22843359b5eSHong Zhang 2295f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(B,maxabsB2,imaxabsB2)); 23043359b5eSHong Zhang 23143359b5eSHong Zhang /* Check maxabsB2 and imaxabsB2 */ 2325f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(maxabsB,&vals)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(maxabsB2,&vals2)); 23443359b5eSHong Zhang for (row=0; row<m; row++) { 23543359b5eSHong Zhang if (PetscAbsScalar(vals[row] - vals2[bs*row]) > PETSC_MACHINE_EPSILON) 23698921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row %" PetscInt_FMT " maxabsB != maxabsB2",row); 23743359b5eSHong Zhang } 2385f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(maxabsB,&vals)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(maxabsB2,&vals2)); 24043359b5eSHong Zhang 24143359b5eSHong Zhang for (col=0; col<n; col++) { 24243359b5eSHong Zhang if (imaxabsB[col] != imaxabsB2[bs*col]/bs) 24398921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"col %" PetscInt_FMT " imaxabsB != imaxabsB2",col); 24443359b5eSHong Zhang } 2455f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&maxabsB)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&maxabsB2)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(imaxabsB,imaxabsB2)); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 2515f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&min)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&max)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&maxabs)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&e)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 256*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 257*b122ec5aSJacob Faibussowitsch return 0; 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown /*TEST 261c4762a1bSJed Brown 262c4762a1bSJed Brown test: 263c4762a1bSJed Brown output_file: output/ex114.out 264c4762a1bSJed Brown 265c4762a1bSJed Brown test: 266c4762a1bSJed Brown suffix: 2 2674e879edeSHong Zhang args: -testcase 1 2684e879edeSHong Zhang output_file: output/ex114.out 269c4762a1bSJed Brown 270c4762a1bSJed Brown test: 271c4762a1bSJed Brown suffix: 3 2724e879edeSHong Zhang args: -testcase 2 2734e879edeSHong Zhang output_file: output/ex114_3.out 2744e879edeSHong Zhang 2754e879edeSHong Zhang test: 2764e879edeSHong Zhang suffix: 4 2774e879edeSHong Zhang args: -testcase 3 2784e879edeSHong Zhang output_file: output/ex114_4.out 2794e879edeSHong Zhang 2804e879edeSHong Zhang test: 2814e879edeSHong Zhang suffix: 5 2824e879edeSHong Zhang nsize: 3 2834e879edeSHong Zhang args: -testcase 0 2844e879edeSHong Zhang output_file: output/ex114_5.out 2854e879edeSHong Zhang 2864e879edeSHong Zhang test: 2874e879edeSHong Zhang suffix: 6 2884e879edeSHong Zhang nsize: 3 2894e879edeSHong Zhang args: -testcase 1 2904e879edeSHong Zhang output_file: output/ex114_6.out 2914e879edeSHong Zhang 2924e879edeSHong Zhang test: 2934e879edeSHong Zhang suffix: 7 2944e879edeSHong Zhang nsize: 3 2954e879edeSHong Zhang args: -testcase 2 2964e879edeSHong Zhang output_file: output/ex114_7.out 2974e879edeSHong Zhang 2984e879edeSHong Zhang test: 2994e879edeSHong Zhang suffix: 8 3004e879edeSHong Zhang nsize: 3 3014e879edeSHong Zhang args: -testcase 3 3024e879edeSHong Zhang output_file: output/ex114_8.out 303c4762a1bSJed Brown 304c4762a1bSJed Brown TEST*/ 305