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 PetscErrorCode ierr; 16c4762a1bSJed Brown IS ip, isrow, iscol; 17c4762a1bSJed Brown PetscRandom rdm; 18c4762a1bSJed Brown PetscBool reorder=PETSC_FALSE; 19c4762a1bSJed Brown MatInfo minfo1,minfo2; 20c4762a1bSJed Brown PetscReal norm1,norm2,tol=10*PETSC_SMALL; 21c4762a1bSJed Brown 22c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 23ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 24*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 25c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);CHKERRQ(ierr); 27c4762a1bSJed Brown 28c4762a1bSJed Brown n = mbs*bs; 29c4762a1bSJed Brown ierr = MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL, &A);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = MatCreateSeqSBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL, &sA);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = MatSetOption(sA,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 35c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&Ii,&J);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = MatGetOwnershipRange(sA,&i,&j);CHKERRQ(ierr); 37c4762a1bSJed Brown if (i-Ii || j-J) { 38c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n");CHKERRQ(ierr); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Assemble matrix */ 42c4762a1bSJed Brown if (bs == 1) { 43c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL);CHKERRQ(ierr); 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; 48c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 54c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 60c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 62c4762a1bSJed Brown } else if (prob ==2) { /* matrix for the five point stencil */ 63c4762a1bSJed Brown n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001); 64*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(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; 70c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown if (i<n1-1) { 74c4762a1bSJed Brown J = Ii + n1; 75c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown if (j>0) { 79c4762a1bSJed Brown J = Ii - 1; 80c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown if (j<n1-1) { 84c4762a1bSJed Brown J = Ii + 1; 85c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 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; 98c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 104c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 110c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 118c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 120c4762a1bSJed Brown col[0]=i; row=i+bs; 121c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } 125c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 127c4762a1bSJed Brown 128c4762a1bSJed Brown ierr = MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* Test MatNorm() */ 132c4762a1bSJed Brown ierr = MatNorm(A,NORM_FROBENIUS,&norm1);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = MatNorm(sA,NORM_FROBENIUS,&norm2);CHKERRQ(ierr); 134c4762a1bSJed Brown norm1 -= norm2; 135c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 1368fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), fnorm1-fnorm2=%16.14e\n",(double)norm1);CHKERRQ(ierr); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown ierr = MatNorm(A,NORM_INFINITY,&norm1);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = MatNorm(sA,NORM_INFINITY,&norm2);CHKERRQ(ierr); 140c4762a1bSJed Brown norm1 -= norm2; 141c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 1428fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), inf_norm1-inf_norm2=%16.14e\n",(double)norm1);CHKERRQ(ierr); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */ 146c4762a1bSJed Brown ierr = MatGetInfo(A,MAT_LOCAL,&minfo1);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = MatGetInfo(sA,MAT_LOCAL,&minfo2);CHKERRQ(ierr); 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) { 151c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatGetInfo()\n");CHKERRQ(ierr); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154c4762a1bSJed Brown ierr = MatGetSize(A,&Ii,&J);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = MatGetSize(sA,&i,&j);CHKERRQ(ierr); 156c4762a1bSJed Brown if (i-Ii || j-J) { 1572f613bf5SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()\n");CHKERRQ(ierr); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown ierr = MatGetBlockSize(A, &Ii);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = MatGetBlockSize(sA, &i);CHKERRQ(ierr); 162c4762a1bSJed Brown if (i-Ii) { 163c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()\n");CHKERRQ(ierr); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */ 167c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&rdm);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,n,&x);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = VecDuplicate(x,&s1);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = VecDuplicate(x,&s2);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = VecDuplicate(x,&b);CHKERRQ(ierr); 174c4762a1bSJed Brown 175c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 176c4762a1bSJed Brown 177c4762a1bSJed Brown ierr = MatDiagonalScale(A,x,x);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = MatDiagonalScale(sA,x,x);CHKERRQ(ierr); 179c4762a1bSJed Brown 180c4762a1bSJed Brown ierr = MatGetDiagonal(A,s1);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = MatGetDiagonal(sA,s2);CHKERRQ(ierr); 182c4762a1bSJed Brown ierr = VecNorm(s1,NORM_1,&norm1);CHKERRQ(ierr); 183c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&norm2);CHKERRQ(ierr); 184c4762a1bSJed Brown norm1 -= norm2; 185c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 186c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal() \n");CHKERRQ(ierr); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown ierr = MatScale(A,alpha);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = MatScale(sA,alpha);CHKERRQ(ierr); 191c4762a1bSJed Brown 192c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 193c4762a1bSJed Brown for (i=0; i<40; i++) { 194c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = MatMult(A,x,s1);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = MatMult(sA,x,s2);CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = VecNorm(s1,NORM_1,&norm1);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&norm2);CHKERRQ(ierr); 199c4762a1bSJed Brown norm1 -= norm2; 200c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 201c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), MatDiagonalScale() or MatScale()\n");CHKERRQ(ierr); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown for (i=0; i<40; i++) { 206c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = VecSetRandom(y,rdm);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = MatMultAdd(sA,x,y,s2);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = VecNorm(s1,NORM_1,&norm1);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&norm2);CHKERRQ(ierr); 212c4762a1bSJed Brown norm1 -= norm2; 213c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 214c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), MatDiagonalScale() or MatScale() \n");CHKERRQ(ierr); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* Test MatReordering() */ 219c4762a1bSJed Brown ierr = MatGetOrdering(A,MATORDERINGNATURAL,&isrow,&iscol);CHKERRQ(ierr); 220c4762a1bSJed Brown ip = isrow; 221c4762a1bSJed Brown 222c4762a1bSJed Brown if (reorder) { 223c4762a1bSJed Brown IS nip; 224c4762a1bSJed Brown PetscInt *nip_ptr; 225c4762a1bSJed Brown ierr = PetscMalloc1(mbs,&nip_ptr);CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = ISGetIndices(ip,&ip_ptr);CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = PetscArraycpy(nip_ptr,ip_ptr,mbs);CHKERRQ(ierr); 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; 230c4762a1bSJed Brown ierr = ISRestoreIndices(ip,&ip_ptr);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,mbs,nip_ptr,PETSC_COPY_VALUES,&nip);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = PetscFree(nip_ptr);CHKERRQ(ierr); 233c4762a1bSJed Brown 234c4762a1bSJed Brown ierr = MatReorderingSeqSBAIJ(sA, ip);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = ISDestroy(&nip);CHKERRQ(ierr); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown ierr = ISDestroy(&iscol);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = ISDestroy(&isrow);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = MatDestroy(&sA);CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = VecDestroy(&s1);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = VecDestroy(&s2);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 248c4762a1bSJed Brown 249c4762a1bSJed Brown ierr = PetscFinalize(); 250c4762a1bSJed Brown return ierr; 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