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