xref: /petsc/src/mat/tests/ex94.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests sequential and parallel MatMatMult() and MatPtAP(), MatTransposeMatMult(), sequential MatMatTransposeMult(), MatRARt()\n\
3c4762a1bSJed Brown Input arguments are:\n\
4c4762a1bSJed Brown   -f0 <input_file> -f1 <input_file> -f2 <input_file> -f3 <input_file> : file to load\n\n";
5c4762a1bSJed Brown /* Example of usage:
6c4762a1bSJed Brown    ./ex94 -f0 <A_binary> -f1 <B_binary> -matmatmult_mat_view ascii::ascii_info -matmatmulttr_mat_view
7c4762a1bSJed Brown    mpiexec -n 3 ./ex94 -f0 medium -f1 medium -f2 arco1 -f3 arco1 -matmatmult_mat_view
8c4762a1bSJed Brown */
9c4762a1bSJed Brown 
10c4762a1bSJed Brown #include <petscmat.h>
11c4762a1bSJed Brown 
12c4762a1bSJed Brown /*
13c4762a1bSJed Brown      B = A - B
14c4762a1bSJed Brown      norm = norm(B)
15c4762a1bSJed Brown */
16c4762a1bSJed Brown PetscErrorCode MatNormDifference(Mat A,Mat B,PetscReal *norm)
17c4762a1bSJed Brown {
18c4762a1bSJed Brown   PetscFunctionBegin;
19*9566063dSJacob Faibussowitsch   PetscCall(MatAXPY(B,-1.0,A,DIFFERENT_NONZERO_PATTERN));
20*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(B,NORM_FROBENIUS,norm));
21c4762a1bSJed Brown   PetscFunctionReturn(0);
22c4762a1bSJed Brown }
23c4762a1bSJed Brown 
24c4762a1bSJed Brown int main(int argc,char **args)
25c4762a1bSJed Brown {
26c4762a1bSJed Brown   Mat            A,A_save,B,AT,ATT,BT,BTT,P,R,C,C1;
27c4762a1bSJed Brown   Vec            x,v1,v2,v3,v4;
28c4762a1bSJed Brown   PetscViewer    viewer;
29c4762a1bSJed Brown   PetscMPIInt    size,rank;
30c4762a1bSJed Brown   PetscInt       i,m,n,j,*idxn,M,N,nzp,rstart,rend;
31c4762a1bSJed Brown   PetscReal      norm,norm_abs,norm_tmp,fill=4.0;
32c4762a1bSJed Brown   PetscRandom    rdm;
33c4762a1bSJed Brown   char           file[4][128];
34c4762a1bSJed Brown   PetscBool      flg,preload = PETSC_TRUE;
35c4762a1bSJed Brown   PetscScalar    *a,rval,alpha,none = -1.0;
36c4762a1bSJed Brown   PetscBool      Test_MatMatMult=PETSC_TRUE,Test_MatMatTr=PETSC_TRUE,Test_MatPtAP=PETSC_TRUE,Test_MatRARt=PETSC_TRUE,Test_MatMatMatMult=PETSC_TRUE;
37c4762a1bSJed Brown   PetscBool      Test_MatAXPY=PETSC_FALSE,view=PETSC_FALSE;
38c4762a1bSJed Brown   PetscInt       pm,pn,pM,pN;
39c4762a1bSJed Brown   MatInfo        info;
40c4762a1bSJed Brown   PetscBool      seqaij;
41c4762a1bSJed Brown   MatType        mattype;
42c4762a1bSJed Brown   Mat            Cdensetest,Pdense,Cdense,Adense;
43c4762a1bSJed Brown 
44*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
45*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
46*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
47c4762a1bSJed Brown 
48*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL));
49*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-matops_view",&view,NULL));
50*9566063dSJacob Faibussowitsch   if (view) PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /*  Load the matrices A_save and B */
53*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f0",file[0],sizeof(file[0]),&flg));
5428b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate a file name for small matrix A with the -f0 option.");
55*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f1",file[1],sizeof(file[1]),&flg));
5628b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate a file name for small matrix B with the -f1 option.");
57*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f2",file[2],sizeof(file[2]),&flg));
58c4762a1bSJed Brown   if (!flg) {
59c4762a1bSJed Brown     preload = PETSC_FALSE;
60c4762a1bSJed Brown   } else {
61*9566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL,NULL,"-f3",file[3],sizeof(file[3]),&flg));
6228b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate a file name for test matrix B with the -f3 option.");
63c4762a1bSJed Brown   }
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   PetscPreLoadBegin(preload,"Load system");
66*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2*PetscPreLoadIt],FILE_MODE_READ,&viewer));
67*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A_save));
68*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A_save));
69*9566063dSJacob Faibussowitsch   PetscCall(MatLoad(A_save,viewer));
70*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
71c4762a1bSJed Brown 
72*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2*PetscPreLoadIt+1],FILE_MODE_READ,&viewer));
73*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
74*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
75*9566063dSJacob Faibussowitsch   PetscCall(MatLoad(B,viewer));
76*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
77c4762a1bSJed Brown 
78*9566063dSJacob Faibussowitsch   PetscCall(MatGetType(B,&mattype));
79c4762a1bSJed Brown 
80*9566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B,&M,&N));
81c4762a1bSJed Brown   nzp  = PetscMax((PetscInt)(0.1*M),5);
82*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc((nzp+1)*(sizeof(PetscInt)+sizeof(PetscScalar)),&idxn));
83c4762a1bSJed Brown   a    = (PetscScalar*)(idxn + nzp);
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* Create vectors v1 and v2 that are compatible with A_save */
86*9566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&v1));
87*9566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A_save,&m,NULL));
88*9566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(v1,m,PETSC_DECIDE));
89*9566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(v1));
90*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(v1,&v2));
91c4762a1bSJed Brown 
92*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm));
93*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
94*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Test MatAXPY()    */
97c4762a1bSJed Brown   /*-------------------*/
98*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-test_MatAXPY",&Test_MatAXPY));
99c4762a1bSJed Brown   if (Test_MatAXPY) {
100c4762a1bSJed Brown     Mat Btmp;
101*9566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A_save,MAT_COPY_VALUES,&A));
102*9566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&Btmp));
103*9566063dSJacob Faibussowitsch     PetscCall(MatAXPY(A,-1.0,B,DIFFERENT_NONZERO_PATTERN)); /* A = -B + A_save */
104c4762a1bSJed Brown 
105*9566063dSJacob Faibussowitsch     PetscCall(MatScale(A,-1.0)); /* A = -A = B - A_save */
106*9566063dSJacob Faibussowitsch     PetscCall(MatAXPY(Btmp,-1.0,A,DIFFERENT_NONZERO_PATTERN)); /* Btmp = -A + B = A_save */
107*9566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A_save,Btmp,10,&flg));
10828b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatAXPY() is incorrect");
109*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
110*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Btmp));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown     Test_MatMatMult    = PETSC_FALSE;
113c4762a1bSJed Brown     Test_MatMatTr      = PETSC_FALSE;
114c4762a1bSJed Brown     Test_MatPtAP       = PETSC_FALSE;
115c4762a1bSJed Brown     Test_MatRARt       = PETSC_FALSE;
116c4762a1bSJed Brown     Test_MatMatMatMult = PETSC_FALSE;
117c4762a1bSJed Brown   }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* 1) Test MatMatMult() */
120c4762a1bSJed Brown   /* ---------------------*/
121c4762a1bSJed Brown   if (Test_MatMatMult) {
122*9566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A_save,MAT_COPY_VALUES,&A));
123*9566063dSJacob Faibussowitsch     PetscCall(MatCreateTranspose(A,&AT));
124*9566063dSJacob Faibussowitsch     PetscCall(MatCreateTranspose(AT,&ATT));
125*9566063dSJacob Faibussowitsch     PetscCall(MatCreateTranspose(B,&BT));
126*9566063dSJacob Faibussowitsch     PetscCall(MatCreateTranspose(BT,&BTT));
127c20d7725SJed Brown 
128*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(AT,B,MAT_INITIAL_MATRIX,fill,&C));
129*9566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(AT,B,C,10,&flg));
13028b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult() for C=AT*B");
131*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
132c20d7725SJed Brown 
133*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(ATT,B,MAT_INITIAL_MATRIX,fill,&C));
134*9566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(ATT,B,C,10,&flg));
13528b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult() for C=ATT*B");
136*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
137c20d7725SJed Brown 
138*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C));
139*9566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(A,B,C,10,&flg));
14028b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult() for reuse C=A*B");
141c20d7725SJed Brown     /* ATT has different matrix type as A (although they have same internal data structure),
142c20d7725SJed Brown        we cannot call MatProductReplaceMats(ATT,NULL,NULL,C) and MatMatMult(ATT,B,MAT_REUSE_MATRIX,fill,&C) */
143*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
144c20d7725SJed Brown 
145*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A,BTT,MAT_INITIAL_MATRIX,fill,&C));
146*9566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(A,BTT,C,10,&flg));
14728b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult() for C=A*BTT");
148*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
149c20d7725SJed Brown 
150*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(ATT,BTT,MAT_INITIAL_MATRIX,fill,&C));
151*9566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(A,B,C,10,&flg));
15228b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult()");
153*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
154c20d7725SJed Brown 
155*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&BTT));
156*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&BT));
157*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&ATT));
158*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&AT));
159c20d7725SJed Brown 
160*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C));
161*9566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(C,"matmatmult_")); /* enable option '-matmatmult_' for matrix C */
162*9566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(C,MAT_GLOBAL_SUM,&info));
163c4762a1bSJed Brown 
164c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX - reuse symbolic C */
165c4762a1bSJed Brown     alpha=1.0;
166c4762a1bSJed Brown     for (i=0; i<2; i++) {
167c4762a1bSJed Brown       alpha -=0.1;
168*9566063dSJacob Faibussowitsch       PetscCall(MatScale(A,alpha));
169*9566063dSJacob Faibussowitsch       PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
170c4762a1bSJed Brown     }
171*9566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(A,B,C,10,&flg));
17228b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult()");
173*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
174c4762a1bSJed Brown 
175c4762a1bSJed Brown     /* Test MatDuplicate() of C=A*B */
176*9566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&C1));
177*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C1));
178*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
179c4762a1bSJed Brown   } /* if (Test_MatMatMult) */
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* 2) Test MatTransposeMatMult() and MatMatTransposeMult() */
182c4762a1bSJed Brown   /* ------------------------------------------------------- */
183c4762a1bSJed Brown   if (Test_MatMatTr) {
184c4762a1bSJed Brown     /* Create P */
185c4762a1bSJed Brown     PetscInt PN,rstart,rend;
186c4762a1bSJed Brown     PN   = M/2;
187c4762a1bSJed Brown     nzp  = 5; /* num of nonzeros in each row of P */
188*9566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&P));
189*9566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(P,PETSC_DECIDE,PETSC_DECIDE,M,PN));
190*9566063dSJacob Faibussowitsch     PetscCall(MatSetType(P,mattype));
191*9566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(P,nzp,NULL));
192*9566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(P,nzp,NULL,nzp,NULL));
193*9566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(P,&rstart,&rend));
194c4762a1bSJed Brown     for (i=0; i<nzp; i++) {
195*9566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm,&a[i]));
196c4762a1bSJed Brown     }
197c4762a1bSJed Brown     for (i=rstart; i<rend; i++) {
198c4762a1bSJed Brown       for (j=0; j<nzp; j++) {
199*9566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValue(rdm,&rval));
200c4762a1bSJed Brown         idxn[j] = (PetscInt)(PetscRealPart(rval)*PN);
201c4762a1bSJed Brown       }
202*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(P,1,&i,nzp,idxn,a,ADD_VALUES));
203c4762a1bSJed Brown     }
204*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY));
205*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY));
206c4762a1bSJed Brown 
207c4762a1bSJed Brown     /* Create R = P^T */
208*9566063dSJacob Faibussowitsch     PetscCall(MatTranspose(P,MAT_INITIAL_MATRIX,&R));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown     { /* Test R = P^T, C1 = R*B */
211*9566063dSJacob Faibussowitsch       PetscCall(MatMatMult(R,B,MAT_INITIAL_MATRIX,fill,&C1));
212*9566063dSJacob Faibussowitsch       PetscCall(MatTranspose(P,MAT_REUSE_MATRIX,&R));
213*9566063dSJacob Faibussowitsch       PetscCall(MatMatMult(R,B,MAT_REUSE_MATRIX,fill,&C1));
214*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&C1));
215c4762a1bSJed Brown     }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown     /* C = P^T*B */
218*9566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(P,B,MAT_INITIAL_MATRIX,fill,&C));
219*9566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(C,MAT_GLOBAL_SUM,&info));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX - reuse symbolic C */
222*9566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(P,B,MAT_REUSE_MATRIX,fill,&C));
223c4762a1bSJed Brown     if (view) {
224*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"C = P^T * B:\n"));
225*9566063dSJacob Faibussowitsch       PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
226c4762a1bSJed Brown     }
227*9566063dSJacob Faibussowitsch     PetscCall(MatProductClear(C));
228c4762a1bSJed Brown     if (view) {
229*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nC = P^T * B after MatProductClear():\n"));
230*9566063dSJacob Faibussowitsch       PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
231c4762a1bSJed Brown     }
232c4762a1bSJed Brown 
233c4762a1bSJed Brown     /* Compare P^T*B and R*B */
234*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(R,B,MAT_INITIAL_MATRIX,fill,&C1));
235*9566063dSJacob Faibussowitsch     PetscCall(MatNormDifference(C,C1,&norm));
2362c71b3e2SJacob Faibussowitsch     PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatTransposeMatMult(): %g",(double)norm);
237*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C1));
238c4762a1bSJed Brown 
239c4762a1bSJed Brown     /* Test MatDuplicate() of C=P^T*B */
240*9566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&C1));
241*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C1));
242*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
243c4762a1bSJed Brown 
244c4762a1bSJed Brown     /* C = B*R^T */
245*9566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&seqaij));
246c4762a1bSJed Brown     if (size == 1 && seqaij) {
247*9566063dSJacob Faibussowitsch       PetscCall(MatMatTransposeMult(B,R,MAT_INITIAL_MATRIX,fill,&C));
248*9566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(C,"matmatmulttr_")); /* enable '-matmatmulttr_' for matrix C */
249*9566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(C,MAT_GLOBAL_SUM,&info));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown       /* Test MAT_REUSE_MATRIX - reuse symbolic C */
252*9566063dSJacob Faibussowitsch       PetscCall(MatMatTransposeMult(B,R,MAT_REUSE_MATRIX,fill,&C));
253c4762a1bSJed Brown 
254c4762a1bSJed Brown       /* Check */
255*9566063dSJacob Faibussowitsch       PetscCall(MatMatMult(B,P,MAT_INITIAL_MATRIX,fill,&C1));
256*9566063dSJacob Faibussowitsch       PetscCall(MatNormDifference(C,C1,&norm));
2572c71b3e2SJacob Faibussowitsch       PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatMatTransposeMult() %g",(double)norm);
258*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&C1));
259*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&C));
260c4762a1bSJed Brown     }
261*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&P));
262*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&R));
263c4762a1bSJed Brown   }
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /* 3) Test MatPtAP() */
266c4762a1bSJed Brown   /*-------------------*/
267c4762a1bSJed Brown   if (Test_MatPtAP) {
268c4762a1bSJed Brown     PetscInt  PN;
269c4762a1bSJed Brown     Mat       Cdup;
270c4762a1bSJed Brown 
271*9566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A_save,MAT_COPY_VALUES,&A));
272*9566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A,&M,&N));
273*9566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A,&m,&n));
274c4762a1bSJed Brown 
275c4762a1bSJed Brown     PN   = M/2;
276c4762a1bSJed Brown     nzp  = (PetscInt)(0.1*PN+1); /* num of nozeros in each row of P */
277*9566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&P));
278*9566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(P,PETSC_DECIDE,PETSC_DECIDE,N,PN));
279*9566063dSJacob Faibussowitsch     PetscCall(MatSetType(P,mattype));
280*9566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(P,nzp,NULL));
281*9566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(P,nzp,NULL,nzp,NULL));
282c4762a1bSJed Brown     for (i=0; i<nzp; i++) {
283*9566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm,&a[i]));
284c4762a1bSJed Brown     }
285*9566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(P,&rstart,&rend));
286c4762a1bSJed Brown     for (i=rstart; i<rend; i++) {
287c4762a1bSJed Brown       for (j=0; j<nzp; j++) {
288*9566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValue(rdm,&rval));
289c4762a1bSJed Brown         idxn[j] = (PetscInt)(PetscRealPart(rval)*PN);
290c4762a1bSJed Brown       }
291*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(P,1,&i,nzp,idxn,a,ADD_VALUES));
292c4762a1bSJed Brown     }
293*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY));
294*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY));
295c4762a1bSJed Brown 
296*9566063dSJacob Faibussowitsch     /* PetscCall(MatView(P,PETSC_VIEWER_STDOUT_WORLD)); */
297*9566063dSJacob Faibussowitsch     PetscCall(MatGetSize(P,&pM,&pN));
298*9566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(P,&pm,&pn));
299*9566063dSJacob Faibussowitsch     PetscCall(MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C));
300c4762a1bSJed Brown 
301c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX - reuse symbolic C */
302c4762a1bSJed Brown     alpha=1.0;
303c4762a1bSJed Brown     for (i=0; i<2; i++) {
304c4762a1bSJed Brown       alpha -=0.1;
305*9566063dSJacob Faibussowitsch       PetscCall(MatScale(A,alpha));
306*9566063dSJacob Faibussowitsch       PetscCall(MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C));
307c4762a1bSJed Brown     }
308c4762a1bSJed Brown 
309c4762a1bSJed Brown     /* Test PtAP ops with P Dense and A either AIJ or SeqDense (it assumes MatPtAP_XAIJ_XAIJ is fine) */
310*9566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&seqaij));
311c4762a1bSJed Brown     if (seqaij) {
312*9566063dSJacob Faibussowitsch       PetscCall(MatConvert(C,MATSEQDENSE,MAT_INITIAL_MATRIX,&Cdensetest));
313*9566063dSJacob Faibussowitsch       PetscCall(MatConvert(P,MATSEQDENSE,MAT_INITIAL_MATRIX,&Pdense));
314c4762a1bSJed Brown     } else {
315*9566063dSJacob Faibussowitsch       PetscCall(MatConvert(C,MATMPIDENSE,MAT_INITIAL_MATRIX,&Cdensetest));
316*9566063dSJacob Faibussowitsch       PetscCall(MatConvert(P,MATMPIDENSE,MAT_INITIAL_MATRIX,&Pdense));
317c4762a1bSJed Brown     }
318c4762a1bSJed Brown 
319c20d7725SJed Brown     /* test with A(AIJ), Pdense -- call MatPtAP_Basic() when np>1 */
320*9566063dSJacob Faibussowitsch     PetscCall(MatPtAP(A,Pdense,MAT_INITIAL_MATRIX,fill,&Cdense));
321*9566063dSJacob Faibussowitsch     PetscCall(MatPtAP(A,Pdense,MAT_REUSE_MATRIX,fill,&Cdense));
322*9566063dSJacob Faibussowitsch     PetscCall(MatPtAPMultEqual(A,Pdense,Cdense,10,&flg));
32328b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatPtAP with A AIJ and P Dense");
324*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Cdense));
325c4762a1bSJed Brown 
326c4762a1bSJed Brown     /* test with A SeqDense */
327c4762a1bSJed Brown     if (seqaij) {
328*9566063dSJacob Faibussowitsch       PetscCall(MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&Adense));
329*9566063dSJacob Faibussowitsch       PetscCall(MatPtAP(Adense,Pdense,MAT_INITIAL_MATRIX,fill,&Cdense));
330*9566063dSJacob Faibussowitsch       PetscCall(MatPtAP(Adense,Pdense,MAT_REUSE_MATRIX,fill,&Cdense));
331*9566063dSJacob Faibussowitsch       PetscCall(MatPtAPMultEqual(Adense,Pdense,Cdense,10,&flg));
33228b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatPtAP with A SeqDense and P SeqDense");
333*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Cdense));
334*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Adense));
335c4762a1bSJed Brown     }
336*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Cdensetest));
337*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Pdense));
338c4762a1bSJed Brown 
339c4762a1bSJed Brown     /* Test MatDuplicate() of C=PtAP and MatView(Cdup,...) */
340*9566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&Cdup));
341c4762a1bSJed Brown     if (view) {
342*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nC = P^T * A * P:\n"));
343*9566063dSJacob Faibussowitsch       PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
344c4762a1bSJed Brown 
345*9566063dSJacob Faibussowitsch       PetscCall(MatProductClear(C));
346*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nC = P^T * A * P after MatProductClear():\n"));
347*9566063dSJacob Faibussowitsch       PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
348c4762a1bSJed Brown 
349*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nCdup:\n"));
350*9566063dSJacob Faibussowitsch       PetscCall(MatView(Cdup,PETSC_VIEWER_STDOUT_WORLD));
351c4762a1bSJed Brown     }
352*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Cdup));
353c4762a1bSJed Brown 
354c4762a1bSJed Brown     if (size>1 || !seqaij) Test_MatRARt = PETSC_FALSE;
355c4762a1bSJed Brown     /* 4) Test MatRARt() */
356c4762a1bSJed Brown     /* ----------------- */
357c4762a1bSJed Brown     if (Test_MatRARt) {
358c20d7725SJed Brown       Mat R, RARt, Rdense, RARtdense;
359*9566063dSJacob Faibussowitsch       PetscCall(MatTranspose(P,MAT_INITIAL_MATRIX,&R));
360c20d7725SJed Brown 
361c20d7725SJed Brown       /* Test MatRARt_Basic(), MatMatMatMult_Basic() */
362*9566063dSJacob Faibussowitsch       PetscCall(MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&Rdense));
363*9566063dSJacob Faibussowitsch       PetscCall(MatRARt(A,Rdense,MAT_INITIAL_MATRIX,2.0,&RARtdense));
364*9566063dSJacob Faibussowitsch       PetscCall(MatRARt(A,Rdense,MAT_REUSE_MATRIX,2.0,&RARtdense));
365c20d7725SJed Brown 
366*9566063dSJacob Faibussowitsch       PetscCall(MatConvert(RARtdense,MATAIJ,MAT_INITIAL_MATRIX,&RARt));
367*9566063dSJacob Faibussowitsch       PetscCall(MatNormDifference(C,RARt,&norm));
3682c71b3e2SJacob Faibussowitsch       PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"|PtAP - RARtdense| = %g",(double)norm);
369*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Rdense));
370*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&RARtdense));
371*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&RARt));
372c20d7725SJed Brown 
373c20d7725SJed Brown       /* Test MatRARt() for aij matrices */
374*9566063dSJacob Faibussowitsch       PetscCall(MatRARt(A,R,MAT_INITIAL_MATRIX,2.0,&RARt));
375*9566063dSJacob Faibussowitsch       PetscCall(MatRARt(A,R,MAT_REUSE_MATRIX,2.0,&RARt));
376*9566063dSJacob Faibussowitsch       PetscCall(MatNormDifference(C,RARt,&norm));
3772c71b3e2SJacob Faibussowitsch       PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"|PtAP - RARt| = %g",(double)norm);
378*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&R));
379*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&RARt));
380c4762a1bSJed Brown     }
381c4762a1bSJed Brown 
382c4762a1bSJed Brown     if (Test_MatMatMatMult && size == 1) {
383c4762a1bSJed Brown       Mat       R, RAP;
384*9566063dSJacob Faibussowitsch       PetscCall(MatTranspose(P,MAT_INITIAL_MATRIX,&R));
385*9566063dSJacob Faibussowitsch       PetscCall(MatMatMatMult(R,A,P,MAT_INITIAL_MATRIX,2.0,&RAP));
386*9566063dSJacob Faibussowitsch       PetscCall(MatMatMatMult(R,A,P,MAT_REUSE_MATRIX,2.0,&RAP));
387*9566063dSJacob Faibussowitsch       PetscCall(MatNormDifference(C,RAP,&norm));
3882c71b3e2SJacob Faibussowitsch       PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PtAP != RAP %g",(double)norm);
389*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&R));
390*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&RAP));
391c4762a1bSJed Brown     }
392c4762a1bSJed Brown 
393c4762a1bSJed Brown     /* Create vector x that is compatible with P */
394*9566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
395*9566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(P,&m,&n));
396*9566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(x,n,PETSC_DECIDE));
397*9566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(x));
398c4762a1bSJed Brown 
399*9566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD,&v3));
400*9566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(v3,n,PETSC_DECIDE));
401*9566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(v3));
402*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(v3,&v4));
403c4762a1bSJed Brown 
404c4762a1bSJed Brown     norm = 0.0;
405c4762a1bSJed Brown     for (i=0; i<10; i++) {
406*9566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(x,rdm));
407*9566063dSJacob Faibussowitsch       PetscCall(MatMult(P,x,v1));
408*9566063dSJacob Faibussowitsch       PetscCall(MatMult(A,v1,v2));  /* v2 = A*P*x */
409c4762a1bSJed Brown 
410*9566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(P,v2,v3)); /* v3 = Pt*A*P*x */
411*9566063dSJacob Faibussowitsch       PetscCall(MatMult(C,x,v4));           /* v3 = C*x   */
412*9566063dSJacob Faibussowitsch       PetscCall(VecNorm(v4,NORM_2,&norm_abs));
413*9566063dSJacob Faibussowitsch       PetscCall(VecAXPY(v4,none,v3));
414*9566063dSJacob Faibussowitsch       PetscCall(VecNorm(v4,NORM_2,&norm_tmp));
415c4762a1bSJed Brown 
416c4762a1bSJed Brown       norm_tmp /= norm_abs;
417c4762a1bSJed Brown       if (norm_tmp > norm) norm = norm_tmp;
418c4762a1bSJed Brown     }
4192c71b3e2SJacob Faibussowitsch     PetscCheckFalse(norm >= PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatPtAP(), |v1 - v2|: %g",(double)norm);
420c4762a1bSJed Brown 
421*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
422*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&P));
423*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
424*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v3));
425*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v4));
426*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
427c4762a1bSJed Brown   }
428c4762a1bSJed Brown 
429c4762a1bSJed Brown   /* Destroy objects */
430*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v1));
431*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v2));
432*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
433*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(idxn));
434c4762a1bSJed Brown 
435*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A_save));
436*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
437c4762a1bSJed Brown 
438c4762a1bSJed Brown   PetscPreLoadEnd();
439*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
440b122ec5aSJacob Faibussowitsch   return 0;
441c4762a1bSJed Brown }
442c4762a1bSJed Brown 
443c4762a1bSJed Brown /*TEST
444c4762a1bSJed Brown 
445c4762a1bSJed Brown    test:
446c4762a1bSJed Brown       suffix: 2_mattransposematmult_matmatmult
447c4762a1bSJed Brown       nsize: 3
448dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
449c20d7725SJed Brown       args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -mattransposematmult_via at*b> ex94_2.tmp 2>&1
450c4762a1bSJed Brown 
451c4762a1bSJed Brown    test:
452c4762a1bSJed Brown       suffix: 2_mattransposematmult_scalable
453c4762a1bSJed Brown       nsize: 3
454dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
455c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -mattransposematmult_via scalable> ex94_2.tmp 2>&1
456c4762a1bSJed Brown       output_file: output/ex94_1.out
457c4762a1bSJed Brown 
458c4762a1bSJed Brown    test:
459c4762a1bSJed Brown       suffix: axpy_mpiaij
460dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
461c4762a1bSJed Brown       nsize: 8
462c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY
463c4762a1bSJed Brown       output_file: output/ex94_1.out
464c4762a1bSJed Brown 
465c4762a1bSJed Brown    test:
466c4762a1bSJed Brown       suffix: axpy_mpibaij
467dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
468c4762a1bSJed Brown       nsize: 8
469c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY -mat_type baij
470c4762a1bSJed Brown       output_file: output/ex94_1.out
471c4762a1bSJed Brown 
472c4762a1bSJed Brown    test:
473c4762a1bSJed Brown       suffix: axpy_mpisbaij
474dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
475c4762a1bSJed Brown       nsize: 8
476c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY -mat_type sbaij
477c4762a1bSJed Brown       output_file: output/ex94_1.out
478c4762a1bSJed Brown 
479c4762a1bSJed Brown    test:
480c4762a1bSJed Brown       suffix: matmatmult
481dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
482c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info
483c4762a1bSJed Brown       output_file: output/ex94_1.out
484c4762a1bSJed Brown 
485c4762a1bSJed Brown    test:
486c4762a1bSJed Brown       suffix: matmatmult_2
487dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
488c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -mat_type mpiaij -viewer_binary_skip_info
489c4762a1bSJed Brown       output_file: output/ex94_1.out
490c4762a1bSJed Brown 
491c4762a1bSJed Brown    test:
492c4762a1bSJed Brown       suffix: matmatmult_scalable
493c4762a1bSJed Brown       nsize: 4
494dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
495c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -matmatmult_via scalable
496c4762a1bSJed Brown       output_file: output/ex94_1.out
497c4762a1bSJed Brown 
498c4762a1bSJed Brown    test:
499c4762a1bSJed Brown       suffix: ptap
500c4762a1bSJed Brown       nsize: 3
501dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
502c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -matptap_via scalable
503c4762a1bSJed Brown       output_file: output/ex94_1.out
504c4762a1bSJed Brown 
505c4762a1bSJed Brown    test:
506c4762a1bSJed Brown       suffix: rap
507c4762a1bSJed Brown       nsize: 3
508dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
509c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium
510c4762a1bSJed Brown       output_file: output/ex94_1.out
511c4762a1bSJed Brown 
512c4762a1bSJed Brown    test:
513c4762a1bSJed Brown       suffix: scalable0
514dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
515c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info
516c4762a1bSJed Brown       output_file: output/ex94_1.out
517c4762a1bSJed Brown 
518c4762a1bSJed Brown    test:
519c4762a1bSJed Brown       suffix: scalable1
520dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
521c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info -matptap_via scalable
522c4762a1bSJed Brown       output_file: output/ex94_1.out
523c4762a1bSJed Brown 
524c4762a1bSJed Brown    test:
525c4762a1bSJed Brown       suffix: view
526c4762a1bSJed Brown       nsize: 2
527dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
528c4762a1bSJed Brown       args: -f0 ${DATAFILESPATH}/matrices/tiny -f1 ${DATAFILESPATH}/matrices/tiny -viewer_binary_skip_info -matops_view
529c4762a1bSJed Brown       output_file: output/ex94_2.out
530c4762a1bSJed Brown 
531c4762a1bSJed Brown TEST*/
532