xref: /petsc/src/mat/tests/ex23.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests the use of interface functions for MATIS matrices.\n\
3*c4762a1bSJed Brown This example tests: MatZeroRows(), MatZeroRowsLocal(), MatView(), MatDuplicate(),\n\
4*c4762a1bSJed Brown MatCopy(), MatCreateSubMatrix(), MatGetLocalSubMatrix(), MatAXPY(), MatShift()\n\
5*c4762a1bSJed Brown MatDiagonalSet(), MatTranspose() and MatPtAP(). It also tests some\n\
6*c4762a1bSJed Brown conversion routines.\n";
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown #include <petscmat.h>
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown PetscErrorCode TestMatZeroRows(Mat,Mat,PetscBool,IS,PetscScalar);
11*c4762a1bSJed Brown PetscErrorCode CheckMat(Mat,Mat,PetscBool,const char*);
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown int main(int argc,char **args)
14*c4762a1bSJed Brown {
15*c4762a1bSJed Brown   Mat                    A,B,A2,B2,T;
16*c4762a1bSJed Brown   Mat                    Aee,Aeo,Aoe,Aoo;
17*c4762a1bSJed Brown   Mat                    *mats;
18*c4762a1bSJed Brown   Vec                    x,y;
19*c4762a1bSJed Brown   MatInfo                info;
20*c4762a1bSJed Brown   ISLocalToGlobalMapping cmap,rmap;
21*c4762a1bSJed Brown   IS                     is,is2,reven,rodd,ceven,codd;
22*c4762a1bSJed Brown   IS                     *rows,*cols;
23*c4762a1bSJed Brown   MatType                lmtype;
24*c4762a1bSJed Brown   PetscScalar            diag = 2.;
25*c4762a1bSJed Brown   PetscInt               n,m,i;
26*c4762a1bSJed Brown   PetscInt               rst,ren,cst,cen,nr,nc;
27*c4762a1bSJed Brown   PetscMPIInt            rank,size;
28*c4762a1bSJed Brown   PetscBool              testT,squaretest,isaij;
29*c4762a1bSJed Brown   PetscBool              permute = PETSC_FALSE;
30*c4762a1bSJed Brown   PetscBool              diffmap = PETSC_TRUE, symmetric = PETSC_FALSE;
31*c4762a1bSJed Brown   PetscErrorCode         ierr;
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
34*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
36*c4762a1bSJed Brown   m = n = 2*size;
37*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-symmetric",&symmetric,NULL);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
40*c4762a1bSJed Brown   if (size > 1 && m < 4) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of rows should be larger or equal 4 for parallel runs");
41*c4762a1bSJed Brown   if (size == 1 && m < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of rows should be larger or equal 2 for uniprocessor runs");
42*c4762a1bSJed Brown   if (n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of cols should be larger or equal 2");
43*c4762a1bSJed Brown   if (symmetric) m = n = PetscMax(m,n);
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown   /* create a MATIS matrix */
46*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = MatSetType(A,MATIS);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
50*c4762a1bSJed Brown   /* This is not the proper setting for MATIS for finite elements, it is just used to test the routines */
51*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,n,0,1,&is);CHKERRQ(ierr);
52*c4762a1bSJed Brown   ierr = ISLocalToGlobalMappingCreateIS(is,&cmap);CHKERRQ(ierr);
53*c4762a1bSJed Brown   ierr = ISDestroy(&is);CHKERRQ(ierr);
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-permmap",&permute,NULL);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-diffmap",&diffmap,NULL);CHKERRQ(ierr);
57*c4762a1bSJed Brown   if (!symmetric && (diffmap || m != n)) {
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,m,permute ? m -1 : 0,permute ? -1 : 1,&is);CHKERRQ(ierr);
60*c4762a1bSJed Brown     ierr = ISLocalToGlobalMappingCreateIS(is,&rmap);CHKERRQ(ierr);
61*c4762a1bSJed Brown     ierr = ISDestroy(&is);CHKERRQ(ierr);
62*c4762a1bSJed Brown     if (m==n && !permute) squaretest = PETSC_TRUE;
63*c4762a1bSJed Brown     else squaretest = PETSC_FALSE;
64*c4762a1bSJed Brown   } else {
65*c4762a1bSJed Brown     ierr = PetscObjectReference((PetscObject)cmap);CHKERRQ(ierr);
66*c4762a1bSJed Brown     rmap = cmap;
67*c4762a1bSJed Brown     squaretest = PETSC_TRUE;
68*c4762a1bSJed Brown   }
69*c4762a1bSJed Brown   ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);
70*c4762a1bSJed Brown   ierr = MatISStoreL2L(A,PETSC_FALSE);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = MatISSetPreallocation(A,3,NULL,0,NULL);CHKERRQ(ierr);
72*c4762a1bSJed Brown   for (i=0; i<m; i++) {
73*c4762a1bSJed Brown     PetscScalar v[3];
74*c4762a1bSJed Brown     PetscInt    cols[3];
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown     cols[0] = (i-1+n)%n;
77*c4762a1bSJed Brown     cols[1] = i%n;
78*c4762a1bSJed Brown     cols[2] = (i+1)%n;
79*c4762a1bSJed Brown     v[0]    = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1);
80*c4762a1bSJed Brown     v[1]    =  2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1);
81*c4762a1bSJed Brown     v[2]    = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1);
82*c4762a1bSJed Brown     ierr    = MatSetValuesLocal(A,1,&i,3,cols,v,ADD_VALUES);CHKERRQ(ierr);
83*c4762a1bSJed Brown   }
84*c4762a1bSJed Brown   if (symmetric) {
85*c4762a1bSJed Brown     ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
86*c4762a1bSJed Brown   }
87*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   ierr = MatISGetLocalMat(A,&B);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = MatGetType(B,&lmtype);CHKERRQ(ierr);
92*c4762a1bSJed Brown 
93*c4762a1bSJed Brown   /* test MatGetInfo */
94*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatGetInfo\n");CHKERRQ(ierr);
95*c4762a1bSJed Brown   if (!PetscGlobalRank) {
96*c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
97*c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
98*c4762a1bSJed Brown   }
99*c4762a1bSJed Brown   ierr = MatGetInfo(A,MAT_LOCAL,&info);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"Process  %2d: %D %D %D %D %D\n",PetscGlobalRank,(PetscInt)info.nz_used,
102*c4762a1bSJed Brown                                             (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = MatGetInfo(A,MAT_GLOBAL_MAX,&info);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"GlobalMax  : %D %D %D %D %D\n",(PetscInt)info.nz_used,
106*c4762a1bSJed Brown                                 (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = MatGetInfo(A,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"GlobalSum  : %D %D %D %D %D\n",(PetscInt)info.nz_used,
109*c4762a1bSJed Brown                                 (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs);CHKERRQ(ierr);
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown   /* test MatView */
112*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatView\n");CHKERRQ(ierr);
113*c4762a1bSJed Brown   ierr = MatView(A,NULL);CHKERRQ(ierr);
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown   /* Create a MPIAIJ matrix, same as A */
116*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = MatSetUp(B);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = MatSetLocalToGlobalMapping(B,rmap,cmap);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocation(B,3,NULL,3,NULL);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = MatMPIBAIJSetPreallocation(B,1,3,NULL,3,NULL);CHKERRQ(ierr);
124*c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE)
125*c4762a1bSJed Brown   ierr = MatHYPRESetPreallocation(B,3,NULL,3,NULL);CHKERRQ(ierr);
126*c4762a1bSJed Brown #endif
127*c4762a1bSJed Brown   ierr = MatISSetPreallocation(B,3,NULL,3,NULL);CHKERRQ(ierr);
128*c4762a1bSJed Brown   for (i=0; i<m; i++) {
129*c4762a1bSJed Brown     PetscScalar v[3];
130*c4762a1bSJed Brown     PetscInt    cols[3];
131*c4762a1bSJed Brown 
132*c4762a1bSJed Brown     cols[0] = (i-1+n)%n;
133*c4762a1bSJed Brown     cols[1] = i%n;
134*c4762a1bSJed Brown     cols[2] = (i+1)%n;
135*c4762a1bSJed Brown     v[0]    = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1);
136*c4762a1bSJed Brown     v[1]    =  2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1);
137*c4762a1bSJed Brown     v[2]    = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1);
138*c4762a1bSJed Brown     ierr    = MatSetValuesLocal(B,1,&i,3,cols,v,ADD_VALUES);CHKERRQ(ierr);
139*c4762a1bSJed Brown   }
140*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
141*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
143*c4762a1bSJed Brown   ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown   /* test CheckMat */
146*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test CheckMat\n");CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = CheckMat(A,B,PETSC_FALSE,"CheckMat");CHKERRQ(ierr);
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown   /* test MatDuplicate and MatAXPY */
150*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatDuplicate and MatAXPY\n");CHKERRQ(ierr);
151*c4762a1bSJed Brown   ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = CheckMat(A,A2,PETSC_FALSE,"MatDuplicate and MatAXPY");CHKERRQ(ierr);
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown   /* test MatConvert */
155*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_IS_XAIJ\n");CHKERRQ(ierr);
156*c4762a1bSJed Brown   ierr = MatConvert(A2,MATAIJ,MAT_INITIAL_MATRIX,&B2);CHKERRQ(ierr);
157*c4762a1bSJed Brown   ierr = CheckMat(B,B2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_INITIAL_MATRIX");CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = MatConvert(A2,MATAIJ,MAT_REUSE_MATRIX,&B2);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = CheckMat(B,B2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_REUSE_MATRIX");CHKERRQ(ierr);
160*c4762a1bSJed Brown   ierr = MatConvert(A2,MATAIJ,MAT_INPLACE_MATRIX,&A2);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = CheckMat(B,A2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_INPLACE_MATRIX");CHKERRQ(ierr);
162*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
163*c4762a1bSJed Brown   ierr = MatDestroy(&B2);CHKERRQ(ierr);
164*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_XAIJ_IS\n");CHKERRQ(ierr);
165*c4762a1bSJed Brown   ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr);
166*c4762a1bSJed Brown   ierr = MatConvert(B2,MATIS,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr);
167*c4762a1bSJed Brown   ierr = CheckMat(A,A2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_INITIAL_MATRIX");CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = MatConvert(B2,MATIS,MAT_REUSE_MATRIX,&A2);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = CheckMat(A,A2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_REUSE_MATRIX");CHKERRQ(ierr);
170*c4762a1bSJed Brown   ierr = MatConvert(B2,MATIS,MAT_INPLACE_MATRIX,&B2);CHKERRQ(ierr);
171*c4762a1bSJed Brown   ierr = CheckMat(A,B2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_INPLACE_MATRIX");CHKERRQ(ierr);
172*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
173*c4762a1bSJed Brown   ierr = MatDestroy(&B2);CHKERRQ(ierr);
174*c4762a1bSJed Brown   ierr = PetscStrcmp(lmtype,MATSEQAIJ,&isaij);CHKERRQ(ierr);
175*c4762a1bSJed Brown   if (size == 1 && isaij) { /* tests special code paths in MatConvert_IS_XAIJ */
176*c4762a1bSJed Brown     PetscInt ri, ci, rr[3] = {0,1,0}, cr[4] = {1,2,0,1}, rk[3] = {0,2,1}, ck[4] = {1,0,3,2};
177*c4762a1bSJed Brown 
178*c4762a1bSJed Brown     for (ri = 0; ri < 2; ri++) {
179*c4762a1bSJed Brown       PetscInt *r;
180*c4762a1bSJed Brown 
181*c4762a1bSJed Brown       r = (PetscInt*)(ri == 0 ? rr : rk);
182*c4762a1bSJed Brown       for (ci = 0; ci < 2; ci++) {
183*c4762a1bSJed Brown         PetscInt *c,rb,cb;
184*c4762a1bSJed Brown 
185*c4762a1bSJed Brown         c = (PetscInt*)(ci == 0 ? cr : ck);
186*c4762a1bSJed Brown         for (rb = 1; rb < 4; rb++) {
187*c4762a1bSJed Brown           ierr = ISCreateBlock(PETSC_COMM_SELF,rb,3,r,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
188*c4762a1bSJed Brown           ierr = ISLocalToGlobalMappingCreateIS(is,&rmap);CHKERRQ(ierr);
189*c4762a1bSJed Brown           ierr = ISDestroy(&is);CHKERRQ(ierr);
190*c4762a1bSJed Brown           for (cb = 1; cb < 4; cb++) {
191*c4762a1bSJed Brown             Mat  T,lT,T2;
192*c4762a1bSJed Brown             char testname[256];
193*c4762a1bSJed Brown 
194*c4762a1bSJed Brown             ierr = PetscSNPrintf(testname,sizeof(testname),"MatConvert_IS_XAIJ special case (%D %D, bs %D %D)",ri,ci,rb,cb);CHKERRQ(ierr);
195*c4762a1bSJed Brown             ierr = PetscPrintf(PETSC_COMM_WORLD,"Test %s\n",testname);CHKERRQ(ierr);
196*c4762a1bSJed Brown 
197*c4762a1bSJed Brown             ierr = ISCreateBlock(PETSC_COMM_SELF,cb,4,c,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
198*c4762a1bSJed Brown             ierr = ISLocalToGlobalMappingCreateIS(is,&cmap);CHKERRQ(ierr);
199*c4762a1bSJed Brown             ierr = ISDestroy(&is);CHKERRQ(ierr);
200*c4762a1bSJed Brown 
201*c4762a1bSJed Brown             ierr = MatCreate(PETSC_COMM_SELF,&T);CHKERRQ(ierr);
202*c4762a1bSJed Brown             ierr = MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,rb*3,cb*4);CHKERRQ(ierr);
203*c4762a1bSJed Brown             ierr = MatSetType(T,MATIS);CHKERRQ(ierr);
204*c4762a1bSJed Brown             ierr = MatSetLocalToGlobalMapping(T,rmap,cmap);CHKERRQ(ierr);
205*c4762a1bSJed Brown             ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
206*c4762a1bSJed Brown             ierr = MatISGetLocalMat(T,&lT);CHKERRQ(ierr);
207*c4762a1bSJed Brown             ierr = MatSetType(lT,MATSEQAIJ);CHKERRQ(ierr);
208*c4762a1bSJed Brown             ierr = MatSeqAIJSetPreallocation(lT,cb*4,NULL);CHKERRQ(ierr);
209*c4762a1bSJed Brown             ierr = MatSetRandom(lT,NULL);CHKERRQ(ierr);
210*c4762a1bSJed Brown             ierr = MatConvert(lT,lmtype,MAT_INPLACE_MATRIX,&lT);CHKERRQ(ierr);
211*c4762a1bSJed Brown             ierr = MatISRestoreLocalMat(T,&lT);CHKERRQ(ierr);
212*c4762a1bSJed Brown             ierr = MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
213*c4762a1bSJed Brown             ierr = MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown             ierr = MatConvert(T,MATAIJ,MAT_INITIAL_MATRIX,&T2);CHKERRQ(ierr);
216*c4762a1bSJed Brown             ierr = CheckMat(T,T2,PETSC_TRUE,"MAT_INITIAL_MATRIX");CHKERRQ(ierr);
217*c4762a1bSJed Brown             ierr = MatConvert(T,MATAIJ,MAT_REUSE_MATRIX,&T2);CHKERRQ(ierr);
218*c4762a1bSJed Brown             ierr = CheckMat(T,T2,PETSC_TRUE,"MAT_REUSE_MATRIX");CHKERRQ(ierr);
219*c4762a1bSJed Brown             ierr = MatDestroy(&T2);CHKERRQ(ierr);
220*c4762a1bSJed Brown             ierr = MatDuplicate(T,MAT_COPY_VALUES,&T2);CHKERRQ(ierr);
221*c4762a1bSJed Brown             ierr = MatConvert(T2,MATAIJ,MAT_INPLACE_MATRIX,&T2);CHKERRQ(ierr);
222*c4762a1bSJed Brown             ierr = CheckMat(T,T2,PETSC_TRUE,"MAT_INPLACE_MATRIX");CHKERRQ(ierr);
223*c4762a1bSJed Brown             ierr = MatDestroy(&T);CHKERRQ(ierr);
224*c4762a1bSJed Brown             ierr = MatDestroy(&T2);CHKERRQ(ierr);
225*c4762a1bSJed Brown           }
226*c4762a1bSJed Brown           ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
227*c4762a1bSJed Brown         }
228*c4762a1bSJed Brown       }
229*c4762a1bSJed Brown     }
230*c4762a1bSJed Brown   }
231*c4762a1bSJed Brown 
232*c4762a1bSJed Brown   /* test MatDiagonalScale */
233*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatDiagonalScale\n");CHKERRQ(ierr);
234*c4762a1bSJed Brown   ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
235*c4762a1bSJed Brown   ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = MatCreateVecs(A,&x,&y);CHKERRQ(ierr);
237*c4762a1bSJed Brown   ierr = VecSetRandom(x,NULL);CHKERRQ(ierr);
238*c4762a1bSJed Brown   if (symmetric) {
239*c4762a1bSJed Brown     ierr = VecCopy(x,y);CHKERRQ(ierr);
240*c4762a1bSJed Brown   } else {
241*c4762a1bSJed Brown     ierr = VecSetRandom(y,NULL);CHKERRQ(ierr);
242*c4762a1bSJed Brown     ierr = VecScale(y,8.);CHKERRQ(ierr);
243*c4762a1bSJed Brown   }
244*c4762a1bSJed Brown   ierr = MatDiagonalScale(A2,y,x);CHKERRQ(ierr);
245*c4762a1bSJed Brown   ierr = MatDiagonalScale(B2,y,x);CHKERRQ(ierr);
246*c4762a1bSJed Brown   ierr = CheckMat(A2,B2,PETSC_FALSE,"MatDiagonalScale");CHKERRQ(ierr);
247*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
248*c4762a1bSJed Brown   ierr = MatDestroy(&B2);CHKERRQ(ierr);
249*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
250*c4762a1bSJed Brown   ierr = VecDestroy(&y);CHKERRQ(ierr);
251*c4762a1bSJed Brown 
252*c4762a1bSJed Brown   /* test MatPtAP (A IS and B AIJ) */
253*c4762a1bSJed Brown   if (isaij && m == n) {
254*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatPtAP\n");CHKERRQ(ierr);
255*c4762a1bSJed Brown     ierr = MatISStoreL2L(A,PETSC_TRUE);CHKERRQ(ierr);
256*c4762a1bSJed Brown     ierr = MatPtAP(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&A2);CHKERRQ(ierr);
257*c4762a1bSJed Brown     ierr = MatPtAP(B,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B2);CHKERRQ(ierr);
258*c4762a1bSJed Brown     ierr = CheckMat(A2,B2,PETSC_FALSE,"MatPtAP MAT_INITIAL_MATRIX");CHKERRQ(ierr);
259*c4762a1bSJed Brown     ierr = MatPtAP(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&A2);CHKERRQ(ierr);
260*c4762a1bSJed Brown     ierr = CheckMat(A2,B2,PETSC_FALSE,"MatPtAP MAT_REUSE_MATRIX");CHKERRQ(ierr);
261*c4762a1bSJed Brown     ierr = MatDestroy(&A2);CHKERRQ(ierr);
262*c4762a1bSJed Brown     ierr = MatDestroy(&B2);CHKERRQ(ierr);
263*c4762a1bSJed Brown   }
264*c4762a1bSJed Brown 
265*c4762a1bSJed Brown   /* test MatGetLocalSubMatrix */
266*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatGetLocalSubMatrix\n");CHKERRQ(ierr);
267*c4762a1bSJed Brown   ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&A2);CHKERRQ(ierr);
268*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,m/2+m%2,0,2,&reven);CHKERRQ(ierr);
269*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,m/2,1,2,&rodd);CHKERRQ(ierr);
270*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,n/2+n%2,0,2,&ceven);CHKERRQ(ierr);
271*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,n/2,1,2,&codd);CHKERRQ(ierr);
272*c4762a1bSJed Brown   ierr = MatGetLocalSubMatrix(A2,reven,ceven,&Aee);CHKERRQ(ierr);
273*c4762a1bSJed Brown   ierr = MatGetLocalSubMatrix(A2,reven,codd,&Aeo);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = MatGetLocalSubMatrix(A2,rodd,ceven,&Aoe);CHKERRQ(ierr);
275*c4762a1bSJed Brown   ierr = MatGetLocalSubMatrix(A2,rodd,codd,&Aoo);CHKERRQ(ierr);
276*c4762a1bSJed Brown   for (i=0; i<m; i++) {
277*c4762a1bSJed Brown     PetscInt    j,je,jo,colse[3], colso[3];
278*c4762a1bSJed Brown     PetscScalar ve[3], vo[3];
279*c4762a1bSJed Brown     PetscScalar v[3];
280*c4762a1bSJed Brown     PetscInt    cols[3];
281*c4762a1bSJed Brown 
282*c4762a1bSJed Brown     cols[0] = (i-1+n)%n;
283*c4762a1bSJed Brown     cols[1] = i%n;
284*c4762a1bSJed Brown     cols[2] = (i+1)%n;
285*c4762a1bSJed Brown     v[0]    = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1);
286*c4762a1bSJed Brown     v[1]    =  2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1);
287*c4762a1bSJed Brown     v[2]    = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1);
288*c4762a1bSJed Brown     for (j=0,je=0,jo=0;j<3;j++) {
289*c4762a1bSJed Brown       if (cols[j]%2) {
290*c4762a1bSJed Brown         vo[jo] = v[j];
291*c4762a1bSJed Brown         colso[jo++] = cols[j]/2;
292*c4762a1bSJed Brown       } else {
293*c4762a1bSJed Brown         ve[je] = v[j];
294*c4762a1bSJed Brown         colse[je++] = cols[j]/2;
295*c4762a1bSJed Brown       }
296*c4762a1bSJed Brown     }
297*c4762a1bSJed Brown     if (i%2) {
298*c4762a1bSJed Brown       PetscInt row = i/2;
299*c4762a1bSJed Brown       ierr = MatSetValuesLocal(Aoe,1,&row,je,colse,ve,ADD_VALUES);CHKERRQ(ierr);
300*c4762a1bSJed Brown       ierr = MatSetValuesLocal(Aoo,1,&row,jo,colso,vo,ADD_VALUES);CHKERRQ(ierr);
301*c4762a1bSJed Brown     } else {
302*c4762a1bSJed Brown       PetscInt row = i/2;
303*c4762a1bSJed Brown       ierr = MatSetValuesLocal(Aee,1,&row,je,colse,ve,ADD_VALUES);CHKERRQ(ierr);
304*c4762a1bSJed Brown       ierr = MatSetValuesLocal(Aeo,1,&row,jo,colso,vo,ADD_VALUES);CHKERRQ(ierr);
305*c4762a1bSJed Brown     }
306*c4762a1bSJed Brown   }
307*c4762a1bSJed Brown   ierr = MatRestoreLocalSubMatrix(A2,reven,ceven,&Aee);CHKERRQ(ierr);
308*c4762a1bSJed Brown   ierr = MatRestoreLocalSubMatrix(A2,reven,codd,&Aeo);CHKERRQ(ierr);
309*c4762a1bSJed Brown   ierr = MatRestoreLocalSubMatrix(A2,rodd,ceven,&Aoe);CHKERRQ(ierr);
310*c4762a1bSJed Brown   ierr = MatRestoreLocalSubMatrix(A2,rodd,codd,&Aoo);CHKERRQ(ierr);
311*c4762a1bSJed Brown   ierr = ISDestroy(&reven);CHKERRQ(ierr);
312*c4762a1bSJed Brown   ierr = ISDestroy(&ceven);CHKERRQ(ierr);
313*c4762a1bSJed Brown   ierr = ISDestroy(&rodd);CHKERRQ(ierr);
314*c4762a1bSJed Brown   ierr = ISDestroy(&codd);CHKERRQ(ierr);
315*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
316*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
317*c4762a1bSJed Brown   ierr = MatAXPY(A2,-1.,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
318*c4762a1bSJed Brown   ierr = CheckMat(A2,NULL,PETSC_FALSE,"MatGetLocalSubMatrix");CHKERRQ(ierr);
319*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
320*c4762a1bSJed Brown 
321*c4762a1bSJed Brown   /* test MatConvert_Nest_IS */
322*c4762a1bSJed Brown   testT = PETSC_FALSE;
323*c4762a1bSJed Brown   ierr  = PetscOptionsGetBool(NULL,NULL,"-test_trans",&testT,NULL);CHKERRQ(ierr);
324*c4762a1bSJed Brown 
325*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_Nest_IS\n");CHKERRQ(ierr);
326*c4762a1bSJed Brown   nr   = 2;
327*c4762a1bSJed Brown   nc   = 2;
328*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nr",&nr,NULL);CHKERRQ(ierr);
329*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nc",&nc,NULL);CHKERRQ(ierr);
330*c4762a1bSJed Brown   if (testT) {
331*c4762a1bSJed Brown     ierr = MatGetOwnershipRange(A,&cst,&cen);CHKERRQ(ierr);
332*c4762a1bSJed Brown     ierr = MatGetOwnershipRangeColumn(A,&rst,&ren);CHKERRQ(ierr);
333*c4762a1bSJed Brown   } else {
334*c4762a1bSJed Brown     ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
335*c4762a1bSJed Brown     ierr = MatGetOwnershipRangeColumn(A,&cst,&cen);CHKERRQ(ierr);
336*c4762a1bSJed Brown   }
337*c4762a1bSJed Brown   ierr = PetscMalloc3(nr,&rows,nc,&cols,2*nr*nc,&mats);CHKERRQ(ierr);
338*c4762a1bSJed Brown   for (i=0;i<nr*nc;i++) {
339*c4762a1bSJed Brown     if (testT) {
340*c4762a1bSJed Brown       ierr = MatCreateTranspose(A,&mats[i]);CHKERRQ(ierr);
341*c4762a1bSJed Brown       ierr = MatTranspose(B,MAT_INITIAL_MATRIX,&mats[i+nr*nc]);CHKERRQ(ierr);
342*c4762a1bSJed Brown     } else {
343*c4762a1bSJed Brown       ierr = MatDuplicate(A,MAT_COPY_VALUES,&mats[i]);CHKERRQ(ierr);
344*c4762a1bSJed Brown       ierr = MatDuplicate(B,MAT_COPY_VALUES,&mats[i+nr*nc]);CHKERRQ(ierr);
345*c4762a1bSJed Brown     }
346*c4762a1bSJed Brown   }
347*c4762a1bSJed Brown   for (i=0;i<nr;i++) {
348*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,ren-rst,i+rst,nr,&rows[i]);CHKERRQ(ierr);
349*c4762a1bSJed Brown   }
350*c4762a1bSJed Brown   for (i=0;i<nc;i++) {
351*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,cen-cst,i+cst,nc,&cols[i]);CHKERRQ(ierr);
352*c4762a1bSJed Brown   }
353*c4762a1bSJed Brown   ierr = MatCreateNest(PETSC_COMM_WORLD,nr,rows,nc,cols,mats,&A2);CHKERRQ(ierr);
354*c4762a1bSJed Brown   ierr = MatCreateNest(PETSC_COMM_WORLD,nr,rows,nc,cols,mats+nr*nc,&B2);CHKERRQ(ierr);
355*c4762a1bSJed Brown   for (i=0;i<nr;i++) {
356*c4762a1bSJed Brown     ierr = ISDestroy(&rows[i]);CHKERRQ(ierr);
357*c4762a1bSJed Brown   }
358*c4762a1bSJed Brown   for (i=0;i<nc;i++) {
359*c4762a1bSJed Brown     ierr = ISDestroy(&cols[i]);CHKERRQ(ierr);
360*c4762a1bSJed Brown   }
361*c4762a1bSJed Brown   for (i=0;i<2*nr*nc;i++) {
362*c4762a1bSJed Brown     ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);
363*c4762a1bSJed Brown   }
364*c4762a1bSJed Brown   ierr = PetscFree3(rows,cols,mats);CHKERRQ(ierr);
365*c4762a1bSJed Brown   ierr = MatConvert(B2,MATAIJ,MAT_INITIAL_MATRIX,&T);CHKERRQ(ierr);
366*c4762a1bSJed Brown   ierr = MatDestroy(&B2);CHKERRQ(ierr);
367*c4762a1bSJed Brown   ierr = MatConvert(A2,MATIS,MAT_INITIAL_MATRIX,&B2);CHKERRQ(ierr);
368*c4762a1bSJed Brown   ierr = CheckMat(B2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_INITIAL_MATRIX");CHKERRQ(ierr);
369*c4762a1bSJed Brown   ierr = MatConvert(A2,MATIS,MAT_REUSE_MATRIX,&B2);CHKERRQ(ierr);
370*c4762a1bSJed Brown   ierr = CheckMat(B2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_REUSE_MATRIX");CHKERRQ(ierr);
371*c4762a1bSJed Brown   ierr = MatDestroy(&B2);CHKERRQ(ierr);
372*c4762a1bSJed Brown   ierr = MatConvert(A2,MATIS,MAT_INPLACE_MATRIX,&A2);CHKERRQ(ierr);
373*c4762a1bSJed Brown   ierr = CheckMat(A2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_INPLACE_MATRIX");CHKERRQ(ierr);
374*c4762a1bSJed Brown   ierr = MatDestroy(&T);CHKERRQ(ierr);
375*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
376*c4762a1bSJed Brown 
377*c4762a1bSJed Brown   /* test MatCreateSubMatrix */
378*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatCreateSubMatrix\n");CHKERRQ(ierr);
379*c4762a1bSJed Brown   if (!rank) {
380*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,1,1,1,&is);CHKERRQ(ierr);
381*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,2,0,1,&is2);CHKERRQ(ierr);
382*c4762a1bSJed Brown   } else if (rank == 1) {
383*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,1,0,1,&is);CHKERRQ(ierr);
384*c4762a1bSJed Brown     if (n > 3) {
385*c4762a1bSJed Brown       ierr = ISCreateStride(PETSC_COMM_WORLD,1,3,1,&is2);CHKERRQ(ierr);
386*c4762a1bSJed Brown     } else {
387*c4762a1bSJed Brown       ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is2);CHKERRQ(ierr);
388*c4762a1bSJed Brown     }
389*c4762a1bSJed Brown   } else if (rank == 2 && n > 4) {
390*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is);CHKERRQ(ierr);
391*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,n-4,4,1,&is2);CHKERRQ(ierr);
392*c4762a1bSJed Brown   } else {
393*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is);CHKERRQ(ierr);
394*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is2);CHKERRQ(ierr);
395*c4762a1bSJed Brown   }
396*c4762a1bSJed Brown   ierr = MatCreateSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr);
397*c4762a1bSJed Brown   ierr = MatCreateSubMatrix(B,is,is,MAT_INITIAL_MATRIX,&B2);CHKERRQ(ierr);
398*c4762a1bSJed Brown   ierr = CheckMat(A2,B2,PETSC_TRUE,"first MatCreateSubMatrix");CHKERRQ(ierr);
399*c4762a1bSJed Brown 
400*c4762a1bSJed Brown   ierr = MatCreateSubMatrix(A,is,is,MAT_REUSE_MATRIX,&A2);CHKERRQ(ierr);
401*c4762a1bSJed Brown   ierr = MatCreateSubMatrix(B,is,is,MAT_REUSE_MATRIX,&B2);CHKERRQ(ierr);
402*c4762a1bSJed Brown   ierr = CheckMat(A2,B2,PETSC_FALSE,"reuse MatCreateSubMatrix");CHKERRQ(ierr);
403*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
404*c4762a1bSJed Brown   ierr = MatDestroy(&B2);CHKERRQ(ierr);
405*c4762a1bSJed Brown 
406*c4762a1bSJed Brown   if (!symmetric) {
407*c4762a1bSJed Brown     ierr = MatCreateSubMatrix(A,is,is2,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr);
408*c4762a1bSJed Brown     ierr = MatCreateSubMatrix(B,is,is2,MAT_INITIAL_MATRIX,&B2);CHKERRQ(ierr);
409*c4762a1bSJed Brown     ierr = MatCreateSubMatrix(A,is,is2,MAT_REUSE_MATRIX,&A2);CHKERRQ(ierr);
410*c4762a1bSJed Brown     ierr = MatCreateSubMatrix(B,is,is2,MAT_REUSE_MATRIX,&B2);CHKERRQ(ierr);
411*c4762a1bSJed Brown     ierr = CheckMat(A2,B2,PETSC_FALSE,"second MatCreateSubMatrix");CHKERRQ(ierr);
412*c4762a1bSJed Brown   }
413*c4762a1bSJed Brown 
414*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
415*c4762a1bSJed Brown   ierr = MatDestroy(&B2);CHKERRQ(ierr);
416*c4762a1bSJed Brown   ierr = ISDestroy(&is);CHKERRQ(ierr);
417*c4762a1bSJed Brown   ierr = ISDestroy(&is2);CHKERRQ(ierr);
418*c4762a1bSJed Brown 
419*c4762a1bSJed Brown   /* Create an IS required by MatZeroRows(): just rank zero provides the rows to be eliminated */
420*c4762a1bSJed Brown   if (size > 1) {
421*c4762a1bSJed Brown     if (!rank) {
422*c4762a1bSJed Brown       PetscInt st,len;
423*c4762a1bSJed Brown 
424*c4762a1bSJed Brown       st   = (m+1)/2;
425*c4762a1bSJed Brown       len  = PetscMin(m/2,PetscMax(m-(m+1)/2-1,0));
426*c4762a1bSJed Brown       ierr = ISCreateStride(PETSC_COMM_WORLD,len,st,1,&is);CHKERRQ(ierr);
427*c4762a1bSJed Brown     } else {
428*c4762a1bSJed Brown       ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is);CHKERRQ(ierr);
429*c4762a1bSJed Brown     }
430*c4762a1bSJed Brown   } else {
431*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_WORLD,1,0,1,&is);CHKERRQ(ierr);
432*c4762a1bSJed Brown   }
433*c4762a1bSJed Brown 
434*c4762a1bSJed Brown   if (squaretest) { /* tests for square matrices only, with same maps for rows and columns */
435*c4762a1bSJed Brown     /* test MatDiagonalSet */
436*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatDiagonalSet\n");CHKERRQ(ierr);
437*c4762a1bSJed Brown     ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
438*c4762a1bSJed Brown     ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr);
439*c4762a1bSJed Brown     ierr = MatCreateVecs(A,NULL,&x);CHKERRQ(ierr);
440*c4762a1bSJed Brown     ierr = VecSetRandom(x,NULL);CHKERRQ(ierr);
441*c4762a1bSJed Brown     ierr = MatDiagonalSet(A2,x,INSERT_VALUES);CHKERRQ(ierr);
442*c4762a1bSJed Brown     ierr = MatDiagonalSet(B2,x,INSERT_VALUES);CHKERRQ(ierr);
443*c4762a1bSJed Brown     ierr = CheckMat(A2,B2,PETSC_FALSE,"MatDiagonalSet");CHKERRQ(ierr);
444*c4762a1bSJed Brown     ierr = VecDestroy(&x);CHKERRQ(ierr);
445*c4762a1bSJed Brown     ierr = MatDestroy(&A2);CHKERRQ(ierr);
446*c4762a1bSJed Brown     ierr = MatDestroy(&B2);CHKERRQ(ierr);
447*c4762a1bSJed Brown 
448*c4762a1bSJed Brown     /* test MatShift (MatShift_IS internally uses MatDiagonalSet_IS with ADD_VALUES) */
449*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatShift\n");CHKERRQ(ierr);
450*c4762a1bSJed Brown     ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
451*c4762a1bSJed Brown     ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr);
452*c4762a1bSJed Brown     ierr = MatShift(A2,2.0);CHKERRQ(ierr);
453*c4762a1bSJed Brown     ierr = MatShift(B2,2.0);CHKERRQ(ierr);
454*c4762a1bSJed Brown     ierr = CheckMat(A2,B2,PETSC_FALSE,"MatShift");CHKERRQ(ierr);
455*c4762a1bSJed Brown     ierr = MatDestroy(&A2);CHKERRQ(ierr);
456*c4762a1bSJed Brown     ierr = MatDestroy(&B2);CHKERRQ(ierr);
457*c4762a1bSJed Brown 
458*c4762a1bSJed Brown     /* nonzero diag value is supported for square matrices only */
459*c4762a1bSJed Brown     ierr = TestMatZeroRows(A,B,PETSC_TRUE,is,diag);CHKERRQ(ierr);
460*c4762a1bSJed Brown   }
461*c4762a1bSJed Brown   ierr = TestMatZeroRows(A,B,squaretest,is,0.0);CHKERRQ(ierr);
462*c4762a1bSJed Brown   ierr = ISDestroy(&is);CHKERRQ(ierr);
463*c4762a1bSJed Brown 
464*c4762a1bSJed Brown   /* test MatTranspose */
465*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatTranspose\n");CHKERRQ(ierr);
466*c4762a1bSJed Brown   ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr);
467*c4762a1bSJed Brown   ierr = MatTranspose(B,MAT_INITIAL_MATRIX,&B2);CHKERRQ(ierr);
468*c4762a1bSJed Brown   ierr = CheckMat(A2,B2,PETSC_FALSE,"initial matrix MatTranspose");CHKERRQ(ierr);
469*c4762a1bSJed Brown 
470*c4762a1bSJed Brown   ierr = MatTranspose(A,MAT_REUSE_MATRIX,&A2);CHKERRQ(ierr);
471*c4762a1bSJed Brown   ierr = CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (not in place) MatTranspose");CHKERRQ(ierr);
472*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
473*c4762a1bSJed Brown 
474*c4762a1bSJed Brown   ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
475*c4762a1bSJed Brown   ierr = MatTranspose(A2,MAT_INPLACE_MATRIX,&A2);CHKERRQ(ierr);
476*c4762a1bSJed Brown   ierr = CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (in place) MatTranspose");CHKERRQ(ierr);
477*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
478*c4762a1bSJed Brown 
479*c4762a1bSJed Brown   ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr);
480*c4762a1bSJed Brown   ierr = CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (different type) MatTranspose");CHKERRQ(ierr);
481*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
482*c4762a1bSJed Brown   ierr = MatDestroy(&B2);CHKERRQ(ierr);
483*c4762a1bSJed Brown 
484*c4762a1bSJed Brown   /* test MatISFixLocalEmpty */
485*c4762a1bSJed Brown   if (isaij) {
486*c4762a1bSJed Brown     PetscInt r[2];
487*c4762a1bSJed Brown 
488*c4762a1bSJed Brown     r[0] = 0;
489*c4762a1bSJed Brown     r[1] = PetscMin(m,n)-1;
490*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatISFixLocalEmpty\n");CHKERRQ(ierr);
491*c4762a1bSJed Brown     ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
492*c4762a1bSJed Brown     ierr = MatISFixLocalEmpty(A2,PETSC_TRUE);CHKERRQ(ierr);
493*c4762a1bSJed Brown     ierr = MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
494*c4762a1bSJed Brown     ierr = MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
495*c4762a1bSJed Brown     ierr = CheckMat(A2,B,PETSC_FALSE,"MatISFixLocalEmpty (null)");CHKERRQ(ierr);
496*c4762a1bSJed Brown 
497*c4762a1bSJed Brown     ierr = MatZeroRows(A2,2,r,0.0,NULL,NULL);CHKERRQ(ierr);
498*c4762a1bSJed Brown     ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr);
499*c4762a1bSJed Brown     ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr);
500*c4762a1bSJed Brown     ierr = MatZeroRows(B2,2,r,0.0,NULL,NULL);CHKERRQ(ierr);
501*c4762a1bSJed Brown     ierr = MatISFixLocalEmpty(A2,PETSC_TRUE);CHKERRQ(ierr);
502*c4762a1bSJed Brown     ierr = MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
503*c4762a1bSJed Brown     ierr = MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
504*c4762a1bSJed Brown     ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr);
505*c4762a1bSJed Brown     ierr = CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (rows)");CHKERRQ(ierr);
506*c4762a1bSJed Brown     ierr = MatDestroy(&A2);CHKERRQ(ierr);
507*c4762a1bSJed Brown 
508*c4762a1bSJed Brown     ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
509*c4762a1bSJed Brown     ierr = MatZeroRows(A2,2,r,0.0,NULL,NULL);CHKERRQ(ierr);
510*c4762a1bSJed Brown     ierr = MatTranspose(A2,MAT_INPLACE_MATRIX,&A2);CHKERRQ(ierr);
511*c4762a1bSJed Brown     ierr = MatTranspose(B2,MAT_INPLACE_MATRIX,&B2);CHKERRQ(ierr);
512*c4762a1bSJed Brown     ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr);
513*c4762a1bSJed Brown     ierr = MatISFixLocalEmpty(A2,PETSC_TRUE);CHKERRQ(ierr);
514*c4762a1bSJed Brown     ierr = MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
515*c4762a1bSJed Brown     ierr = MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
516*c4762a1bSJed Brown     ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr);
517*c4762a1bSJed Brown     ierr = CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (cols)");CHKERRQ(ierr);
518*c4762a1bSJed Brown 
519*c4762a1bSJed Brown     ierr = MatDestroy(&A2);CHKERRQ(ierr);
520*c4762a1bSJed Brown     ierr = MatDestroy(&B2);CHKERRQ(ierr);
521*c4762a1bSJed Brown 
522*c4762a1bSJed Brown     if (squaretest) {
523*c4762a1bSJed Brown       ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
524*c4762a1bSJed Brown       ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr);
525*c4762a1bSJed Brown       ierr = MatZeroRowsColumns(A2,2,r,0.0,NULL,NULL);CHKERRQ(ierr);
526*c4762a1bSJed Brown       ierr = MatZeroRowsColumns(B2,2,r,0.0,NULL,NULL);CHKERRQ(ierr);
527*c4762a1bSJed Brown       ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr);
528*c4762a1bSJed Brown       ierr = MatISFixLocalEmpty(A2,PETSC_TRUE);CHKERRQ(ierr);
529*c4762a1bSJed Brown       ierr = MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
530*c4762a1bSJed Brown       ierr = MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
531*c4762a1bSJed Brown       ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr);
532*c4762a1bSJed Brown       ierr = CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (rows+cols)");CHKERRQ(ierr);
533*c4762a1bSJed Brown       ierr = MatDestroy(&A2);CHKERRQ(ierr);
534*c4762a1bSJed Brown       ierr = MatDestroy(&B2);CHKERRQ(ierr);
535*c4762a1bSJed Brown     }
536*c4762a1bSJed Brown 
537*c4762a1bSJed Brown   }
538*c4762a1bSJed Brown 
539*c4762a1bSJed Brown   /* test MatInvertBlockDiagonal
540*c4762a1bSJed Brown        special cases for block-diagonal matrices */
541*c4762a1bSJed Brown   if (m == n) {
542*c4762a1bSJed Brown     ISLocalToGlobalMapping map;
543*c4762a1bSJed Brown     Mat                    Abd,Bbd;
544*c4762a1bSJed Brown     IS                     is,bis;
545*c4762a1bSJed Brown     const PetscScalar      *isbd,*aijbd;
546*c4762a1bSJed Brown     PetscScalar            *vals;
547*c4762a1bSJed Brown     const PetscInt         *sts,*idxs;
548*c4762a1bSJed Brown     PetscInt               *idxs2,diff,perm,nl,bs,st,en,in;
549*c4762a1bSJed Brown     PetscBool              ok;
550*c4762a1bSJed Brown 
551*c4762a1bSJed Brown     for (diff = 0; diff < 3; diff++) {
552*c4762a1bSJed Brown       for (perm = 0; perm < 3; perm++) {
553*c4762a1bSJed Brown         for (bs = 1; bs < 4; bs++) {
554*c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatInvertBlockDiagonal blockdiag %D %D %D %D\n",n,diff,perm,bs);CHKERRQ(ierr);
555*c4762a1bSJed Brown           ierr = PetscMalloc1(bs*bs,&vals);CHKERRQ(ierr);
556*c4762a1bSJed Brown           ierr = MatGetOwnershipRanges(A,&sts);CHKERRQ(ierr);
557*c4762a1bSJed Brown           switch (diff) {
558*c4762a1bSJed Brown           case 1: /* inverted layout by processes */
559*c4762a1bSJed Brown 	    in = 1;
560*c4762a1bSJed Brown             st = sts[size - rank - 1];
561*c4762a1bSJed Brown             en = sts[size - rank];
562*c4762a1bSJed Brown             nl = en - st;
563*c4762a1bSJed Brown             break;
564*c4762a1bSJed Brown           case 2: /* round-robin layout */
565*c4762a1bSJed Brown             in = size;
566*c4762a1bSJed Brown             st = rank;
567*c4762a1bSJed Brown             nl = n/size;
568*c4762a1bSJed Brown             if (rank < n%size) nl++;
569*c4762a1bSJed Brown             break;
570*c4762a1bSJed Brown           default: /* same layout */
571*c4762a1bSJed Brown             in = 1;
572*c4762a1bSJed Brown             st = sts[rank];
573*c4762a1bSJed Brown             en = sts[rank + 1];
574*c4762a1bSJed Brown             nl = en - st;
575*c4762a1bSJed Brown             break;
576*c4762a1bSJed Brown           }
577*c4762a1bSJed Brown           ierr = ISCreateStride(PETSC_COMM_WORLD,nl,st,in,&is);CHKERRQ(ierr);
578*c4762a1bSJed Brown           ierr = ISGetLocalSize(is,&nl);CHKERRQ(ierr);
579*c4762a1bSJed Brown           ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
580*c4762a1bSJed Brown           ierr = PetscMalloc1(nl,&idxs2);CHKERRQ(ierr);
581*c4762a1bSJed Brown           for (i=0;i<nl;i++) {
582*c4762a1bSJed Brown             switch (perm) { /* invert some of the indices */
583*c4762a1bSJed Brown             case 2:
584*c4762a1bSJed Brown               idxs2[i] = rank%2 ? idxs[i] : idxs[nl-i-1];
585*c4762a1bSJed Brown               break;
586*c4762a1bSJed Brown             case 1:
587*c4762a1bSJed Brown               idxs2[i] = rank%2 ? idxs[nl-i-1] : idxs[i];
588*c4762a1bSJed Brown               break;
589*c4762a1bSJed Brown             default:
590*c4762a1bSJed Brown               idxs2[i] = idxs[i];
591*c4762a1bSJed Brown               break;
592*c4762a1bSJed Brown             }
593*c4762a1bSJed Brown           }
594*c4762a1bSJed Brown           ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
595*c4762a1bSJed Brown           ierr = ISCreateBlock(PETSC_COMM_WORLD,bs,nl,idxs2,PETSC_OWN_POINTER,&bis);CHKERRQ(ierr);
596*c4762a1bSJed Brown           ierr = ISLocalToGlobalMappingCreateIS(bis,&map);CHKERRQ(ierr);
597*c4762a1bSJed Brown           ierr = MatCreateIS(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,bs*n,bs*n,map,NULL,&Abd);CHKERRQ(ierr);
598*c4762a1bSJed Brown           ierr = ISLocalToGlobalMappingDestroy(&map);CHKERRQ(ierr);
599*c4762a1bSJed Brown           ierr = MatISSetPreallocation(Abd,bs,NULL,0,NULL);CHKERRQ(ierr);
600*c4762a1bSJed Brown           for (i=0;i<nl;i++) {
601*c4762a1bSJed Brown             PetscInt b1,b2;
602*c4762a1bSJed Brown 
603*c4762a1bSJed Brown             for (b1=0;b1<bs;b1++) {
604*c4762a1bSJed Brown               for (b2=0;b2<bs;b2++) {
605*c4762a1bSJed Brown                 vals[b1*bs + b2] = i*bs*bs + b1*bs + b2 + 1 + (b1 == b2 ? 1.0 : 0);
606*c4762a1bSJed Brown               }
607*c4762a1bSJed Brown             }
608*c4762a1bSJed Brown             ierr = MatSetValuesBlockedLocal(Abd,1,&i,1,&i,vals,INSERT_VALUES);CHKERRQ(ierr);
609*c4762a1bSJed Brown           }
610*c4762a1bSJed Brown           ierr = MatAssemblyBegin(Abd,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
611*c4762a1bSJed Brown           ierr = MatAssemblyEnd(Abd,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
612*c4762a1bSJed Brown           ierr = MatConvert(Abd,MATAIJ,MAT_INITIAL_MATRIX,&Bbd);CHKERRQ(ierr);
613*c4762a1bSJed Brown           ierr = MatInvertBlockDiagonal(Abd,&isbd);CHKERRQ(ierr);
614*c4762a1bSJed Brown           ierr = MatInvertBlockDiagonal(Bbd,&aijbd);CHKERRQ(ierr);
615*c4762a1bSJed Brown           ierr = MatGetLocalSize(Bbd,&nl,NULL);CHKERRQ(ierr);
616*c4762a1bSJed Brown           ok   = PETSC_TRUE;
617*c4762a1bSJed Brown           for (i=0;i<nl/bs;i++) {
618*c4762a1bSJed Brown             PetscInt b1,b2;
619*c4762a1bSJed Brown 
620*c4762a1bSJed Brown             for (b1=0;b1<bs;b1++) {
621*c4762a1bSJed Brown               for (b2=0;b2<bs;b2++) {
622*c4762a1bSJed Brown                 if (PetscAbsScalar(isbd[i*bs*bs+b1*bs + b2]-aijbd[i*bs*bs+b1*bs + b2]) > PETSC_SMALL) ok = PETSC_FALSE;
623*c4762a1bSJed Brown                 if (!ok) {
624*c4762a1bSJed Brown                   ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] ERROR block %d, entry %d %d: %g %g\n",rank,i,b1,b2,isbd[i*bs*bs+b1*bs + b2],aijbd[i*bs*bs+b1*bs + b2]);CHKERRQ(ierr);
625*c4762a1bSJed Brown                   break;
626*c4762a1bSJed Brown                 }
627*c4762a1bSJed Brown               }
628*c4762a1bSJed Brown               if (!ok) break;
629*c4762a1bSJed Brown             }
630*c4762a1bSJed Brown             if (!ok) break;
631*c4762a1bSJed Brown           }
632*c4762a1bSJed Brown           ierr = MatDestroy(&Abd);CHKERRQ(ierr);
633*c4762a1bSJed Brown           ierr = MatDestroy(&Bbd);CHKERRQ(ierr);
634*c4762a1bSJed Brown           ierr = PetscFree(vals);CHKERRQ(ierr);
635*c4762a1bSJed Brown           ierr = ISDestroy(&is);CHKERRQ(ierr);
636*c4762a1bSJed Brown           ierr = ISDestroy(&bis);CHKERRQ(ierr);
637*c4762a1bSJed Brown         }
638*c4762a1bSJed Brown       }
639*c4762a1bSJed Brown     }
640*c4762a1bSJed Brown   }
641*c4762a1bSJed Brown   /* free testing matrices */
642*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
643*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
644*c4762a1bSJed Brown   ierr = PetscFinalize();
645*c4762a1bSJed Brown   return ierr;
646*c4762a1bSJed Brown }
647*c4762a1bSJed Brown 
648*c4762a1bSJed Brown PetscErrorCode CheckMat(Mat A, Mat B, PetscBool usemult, const char* func)
649*c4762a1bSJed Brown {
650*c4762a1bSJed Brown   Mat            Bcheck;
651*c4762a1bSJed Brown   PetscReal      error;
652*c4762a1bSJed Brown   PetscErrorCode ierr;
653*c4762a1bSJed Brown 
654*c4762a1bSJed Brown   PetscFunctionBeginUser;
655*c4762a1bSJed Brown   if (!usemult && B) {
656*c4762a1bSJed Brown     PetscBool hasnorm;
657*c4762a1bSJed Brown 
658*c4762a1bSJed Brown     ierr = MatHasOperation(B,MATOP_NORM,&hasnorm);CHKERRQ(ierr);
659*c4762a1bSJed Brown     if (!hasnorm) usemult = PETSC_TRUE;
660*c4762a1bSJed Brown   }
661*c4762a1bSJed Brown   if (!usemult) {
662*c4762a1bSJed Brown     if (B) {
663*c4762a1bSJed Brown       MatType Btype;
664*c4762a1bSJed Brown 
665*c4762a1bSJed Brown       ierr = MatGetType(B,&Btype);CHKERRQ(ierr);
666*c4762a1bSJed Brown       ierr = MatConvert(A,Btype,MAT_INITIAL_MATRIX,&Bcheck);CHKERRQ(ierr);
667*c4762a1bSJed Brown     } else {
668*c4762a1bSJed Brown       ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&Bcheck);CHKERRQ(ierr);
669*c4762a1bSJed Brown     }
670*c4762a1bSJed Brown     if (B) { /* if B is present, subtract it */
671*c4762a1bSJed Brown       ierr = MatAXPY(Bcheck,-1.,B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
672*c4762a1bSJed Brown     }
673*c4762a1bSJed Brown     ierr = MatNorm(Bcheck,NORM_INFINITY,&error);CHKERRQ(ierr);
674*c4762a1bSJed Brown     if (error > PETSC_SQRT_MACHINE_EPSILON) {
675*c4762a1bSJed Brown       ISLocalToGlobalMapping rl2g,cl2g;
676*c4762a1bSJed Brown 
677*c4762a1bSJed Brown       ierr = PetscObjectSetName((PetscObject)Bcheck,"Assembled Bcheck");CHKERRQ(ierr);
678*c4762a1bSJed Brown       ierr = MatView(Bcheck,NULL);CHKERRQ(ierr);
679*c4762a1bSJed Brown       if (B) {
680*c4762a1bSJed Brown         ierr = PetscObjectSetName((PetscObject)B,"Assembled AIJ");CHKERRQ(ierr);
681*c4762a1bSJed Brown         ierr = MatView(B,NULL);CHKERRQ(ierr);
682*c4762a1bSJed Brown         ierr = MatDestroy(&Bcheck);CHKERRQ(ierr);
683*c4762a1bSJed Brown         ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&Bcheck);CHKERRQ(ierr);
684*c4762a1bSJed Brown         ierr = PetscObjectSetName((PetscObject)Bcheck,"Assembled IS");CHKERRQ(ierr);
685*c4762a1bSJed Brown         ierr = MatView(Bcheck,NULL);CHKERRQ(ierr);
686*c4762a1bSJed Brown       }
687*c4762a1bSJed Brown       ierr = MatDestroy(&Bcheck);CHKERRQ(ierr);
688*c4762a1bSJed Brown       ierr = PetscObjectSetName((PetscObject)A,"MatIS");CHKERRQ(ierr);
689*c4762a1bSJed Brown       ierr = MatView(A,NULL);CHKERRQ(ierr);
690*c4762a1bSJed Brown       ierr = MatGetLocalToGlobalMapping(A,&rl2g,&cl2g);CHKERRQ(ierr);
691*c4762a1bSJed Brown       ierr = ISLocalToGlobalMappingView(rl2g,NULL);CHKERRQ(ierr);
692*c4762a1bSJed Brown       ierr = ISLocalToGlobalMappingView(cl2g,NULL);CHKERRQ(ierr);
693*c4762a1bSJed Brown       SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ERROR ON %s: %g",func,error);
694*c4762a1bSJed Brown     }
695*c4762a1bSJed Brown     ierr = MatDestroy(&Bcheck);CHKERRQ(ierr);
696*c4762a1bSJed Brown   } else {
697*c4762a1bSJed Brown     PetscBool ok,okt;
698*c4762a1bSJed Brown 
699*c4762a1bSJed Brown     ierr = MatMultEqual(A,B,3,&ok);CHKERRQ(ierr);
700*c4762a1bSJed Brown     ierr = MatMultTransposeEqual(A,B,3,&okt);CHKERRQ(ierr);
701*c4762a1bSJed Brown     if (!ok || !okt) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ERROR ON %s: mult ok ?  %d, multtranspose ok ? %d",func,ok,okt);
702*c4762a1bSJed Brown   }
703*c4762a1bSJed Brown   PetscFunctionReturn(0);
704*c4762a1bSJed Brown }
705*c4762a1bSJed Brown 
706*c4762a1bSJed Brown PetscErrorCode TestMatZeroRows(Mat A, Mat Afull, PetscBool squaretest, IS is, PetscScalar diag)
707*c4762a1bSJed Brown {
708*c4762a1bSJed Brown   Mat                    B,Bcheck,B2 = NULL,lB;
709*c4762a1bSJed Brown   Vec                    x = NULL, b = NULL, b2 = NULL;
710*c4762a1bSJed Brown   ISLocalToGlobalMapping l2gr,l2gc;
711*c4762a1bSJed Brown   PetscReal              error;
712*c4762a1bSJed Brown   char                   diagstr[16];
713*c4762a1bSJed Brown   const PetscInt         *idxs;
714*c4762a1bSJed Brown   PetscInt               rst,ren,i,n,N,d;
715*c4762a1bSJed Brown   PetscMPIInt            rank;
716*c4762a1bSJed Brown   PetscBool              miss,haszerorows;
717*c4762a1bSJed Brown   PetscErrorCode         ierr;
718*c4762a1bSJed Brown 
719*c4762a1bSJed Brown   PetscFunctionBeginUser;
720*c4762a1bSJed Brown   if (diag == 0.) {
721*c4762a1bSJed Brown     ierr = PetscStrcpy(diagstr,"zero");CHKERRQ(ierr);
722*c4762a1bSJed Brown   } else {
723*c4762a1bSJed Brown     ierr = PetscStrcpy(diagstr,"nonzero");CHKERRQ(ierr);
724*c4762a1bSJed Brown   }
725*c4762a1bSJed Brown   ierr = ISView(is,NULL);CHKERRQ(ierr);
726*c4762a1bSJed Brown   ierr = MatGetLocalToGlobalMapping(A,&l2gr,&l2gc);CHKERRQ(ierr);
727*c4762a1bSJed Brown   /* tests MatDuplicate and MatCopy */
728*c4762a1bSJed Brown   if (diag == 0.) {
729*c4762a1bSJed Brown     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
730*c4762a1bSJed Brown   } else {
731*c4762a1bSJed Brown     ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&B);CHKERRQ(ierr);
732*c4762a1bSJed Brown     ierr = MatCopy(A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
733*c4762a1bSJed Brown   }
734*c4762a1bSJed Brown   ierr = MatISGetLocalMat(B,&lB);CHKERRQ(ierr);
735*c4762a1bSJed Brown   ierr = MatHasOperation(lB,MATOP_ZERO_ROWS,&haszerorows);CHKERRQ(ierr);
736*c4762a1bSJed Brown   if (squaretest && haszerorows) {
737*c4762a1bSJed Brown 
738*c4762a1bSJed Brown     ierr = MatCreateVecs(B,&x,&b);CHKERRQ(ierr);
739*c4762a1bSJed Brown     ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr);
740*c4762a1bSJed Brown     ierr = VecSetLocalToGlobalMapping(b,l2gr);CHKERRQ(ierr);
741*c4762a1bSJed Brown     ierr = VecSetLocalToGlobalMapping(x,l2gc);CHKERRQ(ierr);
742*c4762a1bSJed Brown     ierr = VecSetRandom(x,NULL);CHKERRQ(ierr);
743*c4762a1bSJed Brown     ierr = VecSetRandom(b,NULL);CHKERRQ(ierr);
744*c4762a1bSJed Brown     /* mimic b[is] = x[is] */
745*c4762a1bSJed Brown     ierr = VecDuplicate(b,&b2);CHKERRQ(ierr);
746*c4762a1bSJed Brown     ierr = VecSetLocalToGlobalMapping(b2,l2gr);CHKERRQ(ierr);
747*c4762a1bSJed Brown     ierr = VecCopy(b,b2);CHKERRQ(ierr);
748*c4762a1bSJed Brown     ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
749*c4762a1bSJed Brown     ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
750*c4762a1bSJed Brown     ierr = VecGetSize(x,&N);CHKERRQ(ierr);
751*c4762a1bSJed Brown     for (i=0;i<n;i++) {
752*c4762a1bSJed Brown       if (0 <= idxs[i] && idxs[i] < N) {
753*c4762a1bSJed Brown         ierr = VecSetValue(b2,idxs[i],diag,INSERT_VALUES);CHKERRQ(ierr);
754*c4762a1bSJed Brown         ierr = VecSetValue(x,idxs[i],1.,INSERT_VALUES);CHKERRQ(ierr);
755*c4762a1bSJed Brown       }
756*c4762a1bSJed Brown     }
757*c4762a1bSJed Brown     ierr = VecAssemblyBegin(b2);CHKERRQ(ierr);
758*c4762a1bSJed Brown     ierr = VecAssemblyEnd(b2);CHKERRQ(ierr);
759*c4762a1bSJed Brown     ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
760*c4762a1bSJed Brown     ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
761*c4762a1bSJed Brown     ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
762*c4762a1bSJed Brown     /*  test ZeroRows on MATIS */
763*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRows (diag %s)\n",diagstr);CHKERRQ(ierr);
764*c4762a1bSJed Brown     ierr = MatZeroRowsIS(B,is,diag,x,b);CHKERRQ(ierr);
765*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRowsColumns (diag %s)\n",diagstr);CHKERRQ(ierr);
766*c4762a1bSJed Brown     ierr = MatZeroRowsColumnsIS(B2,is,diag,NULL,NULL);CHKERRQ(ierr);
767*c4762a1bSJed Brown   } else if (haszerorows) {
768*c4762a1bSJed Brown     /*  test ZeroRows on MATIS */
769*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRows (diag %s)\n",diagstr);CHKERRQ(ierr);
770*c4762a1bSJed Brown     ierr = MatZeroRowsIS(B,is,diag,NULL,NULL);CHKERRQ(ierr);
771*c4762a1bSJed Brown     b = b2 = x = NULL;
772*c4762a1bSJed Brown   } else {
773*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Skipping MatZeroRows (diag %s)\n",diagstr);CHKERRQ(ierr);
774*c4762a1bSJed Brown     b = b2 = x = NULL;
775*c4762a1bSJed Brown   }
776*c4762a1bSJed Brown 
777*c4762a1bSJed Brown   if (squaretest && haszerorows) {
778*c4762a1bSJed Brown     ierr = VecAXPY(b2,-1.,b);CHKERRQ(ierr);
779*c4762a1bSJed Brown     ierr = VecNorm(b2,NORM_INFINITY,&error);CHKERRQ(ierr);
780*c4762a1bSJed Brown     if (error > PETSC_SQRT_MACHINE_EPSILON) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ERROR IN ZEROROWS ON B %g (diag %s)",error,diagstr);
781*c4762a1bSJed Brown   }
782*c4762a1bSJed Brown 
783*c4762a1bSJed Brown   /* test MatMissingDiagonal */
784*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatMissingDiagonal\n");CHKERRQ(ierr);
785*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
786*c4762a1bSJed Brown   ierr = MatMissingDiagonal(B,&miss,&d);CHKERRQ(ierr);
787*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(B,&rst,&ren);CHKERRQ(ierr);
788*c4762a1bSJed Brown   ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
789*c4762a1bSJed Brown   ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "[%d] [%D,%D) Missing %d, row %D (diag %s)\n",rank,rst,ren,(int)miss,d,diagstr);CHKERRQ(ierr);
790*c4762a1bSJed Brown   ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
791*c4762a1bSJed Brown   ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
792*c4762a1bSJed Brown 
793*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
794*c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
795*c4762a1bSJed Brown   ierr = VecDestroy(&b2);CHKERRQ(ierr);
796*c4762a1bSJed Brown 
797*c4762a1bSJed Brown   /* check the result of ZeroRows with that from MPIAIJ routines
798*c4762a1bSJed Brown      assuming that MatConvert_IS_XAIJ and MatZeroRows_MPIAIJ work fine */
799*c4762a1bSJed Brown   if (haszerorows) {
800*c4762a1bSJed Brown     ierr = MatDuplicate(Afull,MAT_COPY_VALUES,&Bcheck);CHKERRQ(ierr);
801*c4762a1bSJed Brown     ierr = MatSetOption(Bcheck,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
802*c4762a1bSJed Brown     ierr = MatZeroRowsIS(Bcheck,is,diag,NULL,NULL);CHKERRQ(ierr);
803*c4762a1bSJed Brown     ierr = CheckMat(B,Bcheck,PETSC_FALSE,"Zerorows");CHKERRQ(ierr);
804*c4762a1bSJed Brown     ierr = MatDestroy(&Bcheck);CHKERRQ(ierr);
805*c4762a1bSJed Brown   }
806*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
807*c4762a1bSJed Brown 
808*c4762a1bSJed Brown   if (B2) { /* test MatZeroRowsColumns */
809*c4762a1bSJed Brown     ierr = MatDuplicate(Afull,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
810*c4762a1bSJed Brown     ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
811*c4762a1bSJed Brown     ierr = MatZeroRowsColumnsIS(B,is,diag,NULL,NULL);CHKERRQ(ierr);
812*c4762a1bSJed Brown     ierr = CheckMat(B2,B,PETSC_FALSE,"MatZeroRowsColumns");CHKERRQ(ierr);
813*c4762a1bSJed Brown     ierr = MatDestroy(&B);CHKERRQ(ierr);
814*c4762a1bSJed Brown     ierr = MatDestroy(&B2);CHKERRQ(ierr);
815*c4762a1bSJed Brown   }
816*c4762a1bSJed Brown   PetscFunctionReturn(0);
817*c4762a1bSJed Brown }
818*c4762a1bSJed Brown 
819*c4762a1bSJed Brown 
820*c4762a1bSJed Brown /*TEST
821*c4762a1bSJed Brown 
822*c4762a1bSJed Brown    test:
823*c4762a1bSJed Brown       args: -test_trans
824*c4762a1bSJed Brown 
825*c4762a1bSJed Brown    test:
826*c4762a1bSJed Brown       suffix: 2
827*c4762a1bSJed Brown       nsize: 4
828*c4762a1bSJed Brown       args: -matis_convert_local_nest -nr 3 -nc 4
829*c4762a1bSJed Brown 
830*c4762a1bSJed Brown    test:
831*c4762a1bSJed Brown       suffix: 3
832*c4762a1bSJed Brown       nsize: 5
833*c4762a1bSJed Brown       args: -m 11 -n 10 -matis_convert_local_nest -nr 2 -nc 1
834*c4762a1bSJed Brown 
835*c4762a1bSJed Brown    test:
836*c4762a1bSJed Brown       suffix: 4
837*c4762a1bSJed Brown       nsize: 6
838*c4762a1bSJed Brown       args: -m 9 -n 12 -test_trans -nr 2 -nc 7
839*c4762a1bSJed Brown 
840*c4762a1bSJed Brown    test:
841*c4762a1bSJed Brown       suffix: 5
842*c4762a1bSJed Brown       nsize: 6
843*c4762a1bSJed Brown       args: -m 12 -n 12 -test_trans -nr 3 -nc 1
844*c4762a1bSJed Brown 
845*c4762a1bSJed Brown    test:
846*c4762a1bSJed Brown       suffix: 6
847*c4762a1bSJed Brown       args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap
848*c4762a1bSJed Brown 
849*c4762a1bSJed Brown    test:
850*c4762a1bSJed Brown       suffix: 7
851*c4762a1bSJed Brown       args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap
852*c4762a1bSJed Brown 
853*c4762a1bSJed Brown    test:
854*c4762a1bSJed Brown       suffix: 8
855*c4762a1bSJed Brown       args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap
856*c4762a1bSJed Brown 
857*c4762a1bSJed Brown    test:
858*c4762a1bSJed Brown       suffix: 9
859*c4762a1bSJed Brown       nsize: 5
860*c4762a1bSJed Brown       args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap
861*c4762a1bSJed Brown 
862*c4762a1bSJed Brown    test:
863*c4762a1bSJed Brown       suffix: 10
864*c4762a1bSJed Brown       nsize: 5
865*c4762a1bSJed Brown       args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap
866*c4762a1bSJed Brown 
867*c4762a1bSJed Brown    testset:
868*c4762a1bSJed Brown       nsize: 5
869*c4762a1bSJed Brown       args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap
870*c4762a1bSJed Brown       output_file: output/ex23_11.out
871*c4762a1bSJed Brown       test:
872*c4762a1bSJed Brown         suffix: vscat_default
873*c4762a1bSJed Brown         args:
874*c4762a1bSJed Brown       test:
875*c4762a1bSJed Brown         suffix: vscat_sf
876*c4762a1bSJed Brown         args: -vecscatter_type sf
877*c4762a1bSJed Brown 
878*c4762a1bSJed Brown    test:
879*c4762a1bSJed Brown       suffix: 12
880*c4762a1bSJed Brown       nsize: 3
881*c4762a1bSJed Brown       args: -m 12 -n 12 -symmetric -matis_localmat_type sbaij -test_trans -nr 2 -nc 3
882*c4762a1bSJed Brown 
883*c4762a1bSJed Brown    testset:
884*c4762a1bSJed Brown       output_file: output/ex23_13.out
885*c4762a1bSJed Brown       nsize: 3
886*c4762a1bSJed Brown       args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -diffmap -permmap
887*c4762a1bSJed Brown       filter: grep -v "type:"
888*c4762a1bSJed Brown       test:
889*c4762a1bSJed Brown         suffix: baij
890*c4762a1bSJed Brown         args: -matis_localmat_type baij
891*c4762a1bSJed Brown       test:
892*c4762a1bSJed Brown         requires: viennacl
893*c4762a1bSJed Brown         suffix: viennacl
894*c4762a1bSJed Brown         args: -matis_localmat_type aijviennacl
895*c4762a1bSJed Brown       test:
896*c4762a1bSJed Brown         requires: cuda
897*c4762a1bSJed Brown         suffix: cusparse
898*c4762a1bSJed Brown         args: -matis_localmat_type aijcusparse
899*c4762a1bSJed Brown       test:
900*c4762a1bSJed Brown         requires: cuda
901*c4762a1bSJed Brown         suffix: cusparse_sf
902*c4762a1bSJed Brown         args: -matis_localmat_type aijcusparse -vecscatter_type sf
903*c4762a1bSJed Brown 
904*c4762a1bSJed Brown TEST*/
905