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