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