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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 22*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 232c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 24*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 25*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown n = mbs*bs; 28*9566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL, &A)); 29*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 30*9566063dSJacob Faibussowitsch PetscCall(MatCreateSeqSBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL, &sA)); 31*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(sA,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 34*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Ii,&J)); 35*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(sA,&i,&j)); 36c4762a1bSJed Brown if (i-Ii || j-J) { 37*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n")); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Assemble matrix */ 41c4762a1bSJed Brown if (bs == 1) { 42*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL)); 43c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */ 44c4762a1bSJed Brown value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 45c4762a1bSJed Brown for (i=1; i<n-1; i++) { 46c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 47*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 48*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1; 51c4762a1bSJed Brown 52c4762a1bSJed Brown value[0]= 0.1; value[1]=-1; value[2]=2; 53*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 54*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 55c4762a1bSJed Brown 56c4762a1bSJed Brown i = 0; col[0] = 0; col[1] = 1; col[2]=n-1; 57c4762a1bSJed Brown 58c4762a1bSJed Brown value[0] = 2.0; value[1] = -1.0; value[2]=0.1; 59*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 60*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 61c4762a1bSJed Brown } else if (prob ==2) { /* matrix for the five point stencil */ 62c4762a1bSJed Brown n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001); 632c71b3e2SJacob Faibussowitsch PetscCheckFalse(n1*n1 - n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive integer!"); 64c4762a1bSJed Brown for (i=0; i<n1; i++) { 65c4762a1bSJed Brown for (j=0; j<n1; j++) { 66c4762a1bSJed Brown Ii = j + n1*i; 67c4762a1bSJed Brown if (i>0) { 68c4762a1bSJed Brown J = Ii - n1; 69*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 70*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown if (i<n1-1) { 73c4762a1bSJed Brown J = Ii + n1; 74*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 75*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown if (j>0) { 78c4762a1bSJed Brown J = Ii - 1; 79*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 80*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown if (j<n1-1) { 83c4762a1bSJed Brown J = Ii + 1; 84*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 85*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown } 88c4762a1bSJed Brown } 89c4762a1bSJed Brown } 90c4762a1bSJed Brown } else { /* bs > 1 */ 91c4762a1bSJed Brown #if defined(DIAGB) 92c4762a1bSJed Brown for (block=0; block<n/bs; block++) { 93c4762a1bSJed Brown /* diagonal blocks */ 94c4762a1bSJed Brown value[0] = -1.0; value[1] = 4.0; value[2] = -1.0; 95c4762a1bSJed Brown for (i=1+block*bs; i<bs-1+block*bs; i++) { 96c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 97*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 98*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs; 101c4762a1bSJed Brown 102c4762a1bSJed Brown value[0]=-1.0; value[1]=4.0; 103*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 104*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs; 107c4762a1bSJed Brown 108c4762a1bSJed Brown value[0]=4.0; value[1] = -1.0; 109*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 110*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES)); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown #endif 113c4762a1bSJed Brown /* off-diagonal blocks */ 114c4762a1bSJed Brown value[0]=-1.0; 115c4762a1bSJed Brown for (i=0; i<(n/bs-1)*bs; i++) { 116c4762a1bSJed Brown col[0]=i+bs; 117*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES)); 118*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES)); 119c4762a1bSJed Brown col[0]=i; row=i+bs; 120*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES)); 121*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES)); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown } 124*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 125*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 126c4762a1bSJed Brown 127*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY)); 128*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* Test MatNorm() */ 131*9566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_FROBENIUS,&norm1)); 132*9566063dSJacob Faibussowitsch PetscCall(MatNorm(sA,NORM_FROBENIUS,&norm2)); 133c4762a1bSJed Brown norm1 -= norm2; 134c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 135*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), fnorm1-fnorm2=%16.14e\n",(double)norm1)); 136c4762a1bSJed Brown } 137*9566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_INFINITY,&norm1)); 138*9566063dSJacob Faibussowitsch PetscCall(MatNorm(sA,NORM_INFINITY,&norm2)); 139c4762a1bSJed Brown norm1 -= norm2; 140c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 141*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), inf_norm1-inf_norm2=%16.14e\n",(double)norm1)); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */ 145*9566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A,MAT_LOCAL,&minfo1)); 146*9566063dSJacob Faibussowitsch PetscCall(MatGetInfo(sA,MAT_LOCAL,&minfo2)); 147c4762a1bSJed Brown i = (int) (minfo1.nz_used - minfo2.nz_used); 148c4762a1bSJed Brown j = (int) (minfo1.nz_allocated - minfo2.nz_allocated); 149c4762a1bSJed Brown if (i<0 || j<0) { 150*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetInfo()\n")); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&Ii,&J)); 154*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(sA,&i,&j)); 155c4762a1bSJed Brown if (i-Ii || j-J) { 156*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()\n")); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159*9566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &Ii)); 160*9566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(sA, &i)); 161c4762a1bSJed Brown if (i-Ii) { 162*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()\n")); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */ 166*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 167*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 168*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&x)); 169*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&s1)); 170*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&s2)); 171*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&y)); 172*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&b)); 173c4762a1bSJed Brown 174*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 175c4762a1bSJed Brown 176*9566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,x,x)); 177*9566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(sA,x,x)); 178c4762a1bSJed Brown 179*9566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A,s1)); 180*9566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(sA,s2)); 181*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_1,&norm1)); 182*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_1,&norm2)); 183c4762a1bSJed Brown norm1 -= norm2; 184c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 185*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal() \n")); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188*9566063dSJacob Faibussowitsch PetscCall(MatScale(A,alpha)); 189*9566063dSJacob Faibussowitsch PetscCall(MatScale(sA,alpha)); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 192c4762a1bSJed Brown for (i=0; i<40; i++) { 193*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 194*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,s1)); 195*9566063dSJacob Faibussowitsch PetscCall(MatMult(sA,x,s2)); 196*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_1,&norm1)); 197*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_1,&norm2)); 198c4762a1bSJed Brown norm1 -= norm2; 199c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 200*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), MatDiagonalScale() or MatScale()\n")); 201c4762a1bSJed Brown } 202c4762a1bSJed Brown } 203c4762a1bSJed Brown 204c4762a1bSJed Brown for (i=0; i<40; i++) { 205*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 206*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y,rdm)); 207*9566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,x,y,s1)); 208*9566063dSJacob Faibussowitsch PetscCall(MatMultAdd(sA,x,y,s2)); 209*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_1,&norm1)); 210*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_1,&norm2)); 211c4762a1bSJed Brown norm1 -= norm2; 212c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 213*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), MatDiagonalScale() or MatScale() \n")); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* Test MatReordering() */ 218*9566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A,MATORDERINGNATURAL,&isrow,&iscol)); 219c4762a1bSJed Brown ip = isrow; 220c4762a1bSJed Brown 221c4762a1bSJed Brown if (reorder) { 222c4762a1bSJed Brown IS nip; 223c4762a1bSJed Brown PetscInt *nip_ptr; 224*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mbs,&nip_ptr)); 225*9566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ip,&ip_ptr)); 226*9566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(nip_ptr,ip_ptr,mbs)); 227c4762a1bSJed Brown i = nip_ptr[1]; nip_ptr[1] = nip_ptr[mbs-2]; nip_ptr[mbs-2] = i; 228c4762a1bSJed Brown i = nip_ptr[0]; nip_ptr[0] = nip_ptr[mbs-1]; nip_ptr[mbs-1] = i; 229*9566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(ip,&ip_ptr)); 230*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,mbs,nip_ptr,PETSC_COPY_VALUES,&nip)); 231*9566063dSJacob Faibussowitsch PetscCall(PetscFree(nip_ptr)); 232c4762a1bSJed Brown 233*9566063dSJacob Faibussowitsch PetscCall(MatReorderingSeqSBAIJ(sA, ip)); 234*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&nip)); 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 238*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 239*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 240*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 241*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 242*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 243*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 244*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 245*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 246*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 247c4762a1bSJed Brown 248*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 249b122ec5aSJacob Faibussowitsch return 0; 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252c4762a1bSJed Brown /*TEST 253c4762a1bSJed Brown 254c4762a1bSJed Brown test: 255c4762a1bSJed Brown args: -bs {{1 2 3 4 5 6 7 8}} 256c4762a1bSJed Brown 257c4762a1bSJed Brown TEST*/ 258