1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests the various sequential routines in MatSBAIJ format. Same as ex74.c except diagonal entries of the matrices are zeros.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Vec x,y,b,s1,s2; 9c4762a1bSJed Brown Mat A; /* linear system matrix */ 10c4762a1bSJed Brown Mat sA; /* symmetric part of the matrices */ 11c4762a1bSJed Brown PetscInt n,mbs=16,bs=1,nz=3,prob=2,i,j,col[3],row,Ii,J,n1; 12c4762a1bSJed Brown const PetscInt *ip_ptr; 13c4762a1bSJed Brown PetscScalar neg_one = -1.0,value[3],alpha=0.1; 14c4762a1bSJed Brown PetscMPIInt size; 15c4762a1bSJed Brown IS ip, isrow, iscol; 16c4762a1bSJed Brown PetscRandom rdm; 17c4762a1bSJed Brown PetscBool reorder=PETSC_FALSE; 18c4762a1bSJed Brown MatInfo minfo1,minfo2; 19c4762a1bSJed Brown PetscReal norm1,norm2,tol=10*PETSC_SMALL; 20c4762a1bSJed Brown 21*327415f7SBarry Smith PetscFunctionBeginUser; 229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 24be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown n = mbs*bs; 299566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL, &A)); 309566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 319566063dSJacob Faibussowitsch PetscCall(MatCreateSeqSBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL, &sA)); 329566063dSJacob Faibussowitsch PetscCall(MatSetOption(sA,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Ii,&J)); 369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(sA,&i,&j)); 37c4762a1bSJed Brown if (i-Ii || j-J) { 389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n")); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Assemble matrix */ 42c4762a1bSJed Brown if (bs == 1) { 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL)); 44c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */ 45c4762a1bSJed Brown value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 46c4762a1bSJed Brown for (i=1; i<n-1; i++) { 47c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 489566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 499566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1; 52c4762a1bSJed Brown 53c4762a1bSJed Brown value[0]= 0.1; value[1]=-1; value[2]=2; 549566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 559566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown i = 0; col[0] = 0; col[1] = 1; col[2]=n-1; 58c4762a1bSJed Brown 59c4762a1bSJed Brown value[0] = 2.0; value[1] = -1.0; value[2]=0.1; 609566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 619566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 62c4762a1bSJed Brown } else if (prob ==2) { /* matrix for the five point stencil */ 63c4762a1bSJed Brown n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001); 64bc30f867SBarry Smith PetscCheck(n1*n1 == n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive integer!"); 65c4762a1bSJed Brown for (i=0; i<n1; i++) { 66c4762a1bSJed Brown for (j=0; j<n1; j++) { 67c4762a1bSJed Brown Ii = j + n1*i; 68c4762a1bSJed Brown if (i>0) { 69c4762a1bSJed Brown J = Ii - n1; 709566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 719566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown if (i<n1-1) { 74c4762a1bSJed Brown J = Ii + n1; 759566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 769566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown if (j>0) { 79c4762a1bSJed Brown J = Ii - 1; 809566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 819566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown if (j<n1-1) { 84c4762a1bSJed Brown J = Ii + 1; 859566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 869566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown } 89c4762a1bSJed Brown } 90c4762a1bSJed Brown } 91c4762a1bSJed Brown } else { /* bs > 1 */ 92c4762a1bSJed Brown #if defined(DIAGB) 93c4762a1bSJed Brown for (block=0; block<n/bs; block++) { 94c4762a1bSJed Brown /* diagonal blocks */ 95c4762a1bSJed Brown value[0] = -1.0; value[1] = 4.0; value[2] = -1.0; 96c4762a1bSJed Brown for (i=1+block*bs; i<bs-1+block*bs; i++) { 97c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 989566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 999566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs; 102c4762a1bSJed Brown 103c4762a1bSJed Brown value[0]=-1.0; value[1]=4.0; 1049566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 1059566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs; 108c4762a1bSJed Brown 109c4762a1bSJed Brown value[0]=4.0; value[1] = -1.0; 1109566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 1119566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES)); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown #endif 114c4762a1bSJed Brown /* off-diagonal blocks */ 115c4762a1bSJed Brown value[0]=-1.0; 116c4762a1bSJed Brown for (i=0; i<(n/bs-1)*bs; i++) { 117c4762a1bSJed Brown col[0]=i+bs; 1189566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES)); 1199566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES)); 120c4762a1bSJed Brown col[0]=i; row=i+bs; 1219566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES)); 1229566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES)); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } 1259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY)); 1299566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY)); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* Test MatNorm() */ 1329566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_FROBENIUS,&norm1)); 1339566063dSJacob Faibussowitsch PetscCall(MatNorm(sA,NORM_FROBENIUS,&norm2)); 134c4762a1bSJed Brown norm1 -= norm2; 135c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 1369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), fnorm1-fnorm2=%16.14e\n",(double)norm1)); 137c4762a1bSJed Brown } 1389566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_INFINITY,&norm1)); 1399566063dSJacob Faibussowitsch PetscCall(MatNorm(sA,NORM_INFINITY,&norm2)); 140c4762a1bSJed Brown norm1 -= norm2; 141c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 1429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), inf_norm1-inf_norm2=%16.14e\n",(double)norm1)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */ 1469566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A,MAT_LOCAL,&minfo1)); 1479566063dSJacob Faibussowitsch PetscCall(MatGetInfo(sA,MAT_LOCAL,&minfo2)); 148c4762a1bSJed Brown i = (int) (minfo1.nz_used - minfo2.nz_used); 149c4762a1bSJed Brown j = (int) (minfo1.nz_allocated - minfo2.nz_allocated); 150c4762a1bSJed Brown if (i<0 || j<0) { 1519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetInfo()\n")); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 1549566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&Ii,&J)); 1559566063dSJacob Faibussowitsch PetscCall(MatGetSize(sA,&i,&j)); 156c4762a1bSJed Brown if (i-Ii || j-J) { 1579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()\n")); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 1609566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &Ii)); 1619566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(sA, &i)); 162c4762a1bSJed Brown if (i-Ii) { 1639566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()\n")); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */ 1679566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 1689566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 1699566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&x)); 1709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&s1)); 1719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&s2)); 1729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&y)); 1739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&b)); 174c4762a1bSJed Brown 1759566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 176c4762a1bSJed Brown 1779566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,x,x)); 1789566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(sA,x,x)); 179c4762a1bSJed Brown 1809566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A,s1)); 1819566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(sA,s2)); 1829566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_1,&norm1)); 1839566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_1,&norm2)); 184c4762a1bSJed Brown norm1 -= norm2; 185c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 1869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal() \n")); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 1899566063dSJacob Faibussowitsch PetscCall(MatScale(A,alpha)); 1909566063dSJacob Faibussowitsch PetscCall(MatScale(sA,alpha)); 191c4762a1bSJed Brown 192c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 193c4762a1bSJed Brown for (i=0; i<40; i++) { 1949566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 1959566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,s1)); 1969566063dSJacob Faibussowitsch PetscCall(MatMult(sA,x,s2)); 1979566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_1,&norm1)); 1989566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_1,&norm2)); 199c4762a1bSJed Brown norm1 -= norm2; 200c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), MatDiagonalScale() or MatScale()\n")); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown for (i=0; i<40; i++) { 2069566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 2079566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y,rdm)); 2089566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,x,y,s1)); 2099566063dSJacob Faibussowitsch PetscCall(MatMultAdd(sA,x,y,s2)); 2109566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_1,&norm1)); 2119566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_1,&norm2)); 212c4762a1bSJed Brown norm1 -= norm2; 213c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), MatDiagonalScale() or MatScale() \n")); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* Test MatReordering() */ 2199566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A,MATORDERINGNATURAL,&isrow,&iscol)); 220c4762a1bSJed Brown ip = isrow; 221c4762a1bSJed Brown 222c4762a1bSJed Brown if (reorder) { 223c4762a1bSJed Brown IS nip; 224c4762a1bSJed Brown PetscInt *nip_ptr; 2259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mbs,&nip_ptr)); 2269566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ip,&ip_ptr)); 2279566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(nip_ptr,ip_ptr,mbs)); 228c4762a1bSJed Brown i = nip_ptr[1]; nip_ptr[1] = nip_ptr[mbs-2]; nip_ptr[mbs-2] = i; 229c4762a1bSJed Brown i = nip_ptr[0]; nip_ptr[0] = nip_ptr[mbs-1]; nip_ptr[mbs-1] = i; 2309566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(ip,&ip_ptr)); 2319566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,mbs,nip_ptr,PETSC_COPY_VALUES,&nip)); 2329566063dSJacob Faibussowitsch PetscCall(PetscFree(nip_ptr)); 233c4762a1bSJed Brown 2349566063dSJacob Faibussowitsch PetscCall(MatReorderingSeqSBAIJ(sA, ip)); 2359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&nip)); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 2389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 2399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 2409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 2429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 2459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 2469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2479566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 248c4762a1bSJed Brown 2499566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 250b122ec5aSJacob Faibussowitsch return 0; 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253c4762a1bSJed Brown /*TEST 254c4762a1bSJed Brown 255c4762a1bSJed Brown test: 256c4762a1bSJed Brown args: -bs {{1 2 3 4 5 6 7 8}} 257c4762a1bSJed Brown 258c4762a1bSJed Brown TEST*/ 259