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