xref: /petsc/src/mat/tests/ex70.c (revision 75ab9b9f8a149e82ac5d3b97f5fa96f50069e065)
16280154eSStefano Zampini #include <petscmat.h>
26280154eSStefano Zampini 
3*75ab9b9fSStefano Zampini static char help[] = "Tests MatMat operations with MAT_REUSE_MATRIX and already allocated dense result.\n\n";
4*75ab9b9fSStefano Zampini 
5*75ab9b9fSStefano Zampini static PetscScalar MAGIC_NUMBER = 12345;
66280154eSStefano Zampini 
76280154eSStefano Zampini static PetscErrorCode CheckLocal(Mat A, Mat B, PetscScalar *a, PetscScalar *b)
86280154eSStefano Zampini {
96280154eSStefano Zampini   PetscErrorCode ierr;
106280154eSStefano Zampini   PetscBool      wA = PETSC_FALSE, wB = PETSC_FALSE;
11*75ab9b9fSStefano Zampini   PetscBool      wAv = PETSC_FALSE, wBv = PETSC_FALSE;
12*75ab9b9fSStefano Zampini   PetscInt       lda,i,j,m,n;
136280154eSStefano Zampini 
146280154eSStefano Zampini   PetscFunctionBegin;
156280154eSStefano Zampini   if (a) {
166280154eSStefano Zampini     const PetscScalar *Aa;
176280154eSStefano Zampini     ierr = MatDenseGetArrayRead(A,&Aa);CHKERRQ(ierr);
186280154eSStefano Zampini     wA   = (PetscBool)(a != Aa);
19*75ab9b9fSStefano Zampini     ierr = MatDenseGetLDA(A,&lda);CHKERRQ(ierr);
20*75ab9b9fSStefano Zampini     ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
21*75ab9b9fSStefano Zampini     for (j=0;j<n;j++) {
22*75ab9b9fSStefano Zampini       for (i=m;i<lda;i++) {
23*75ab9b9fSStefano Zampini         if (Aa[j*lda +i] != MAGIC_NUMBER) wAv = PETSC_TRUE;
24*75ab9b9fSStefano Zampini       }
25*75ab9b9fSStefano Zampini     }
266280154eSStefano Zampini     ierr = MatDenseRestoreArrayRead(A,&Aa);CHKERRQ(ierr);
276280154eSStefano Zampini   }
286280154eSStefano Zampini   if (b) {
296280154eSStefano Zampini     const PetscScalar *Bb;
306280154eSStefano Zampini     ierr = MatDenseGetArrayRead(B,&Bb);CHKERRQ(ierr);
316280154eSStefano Zampini     wB   = (PetscBool)(b != Bb);
32*75ab9b9fSStefano Zampini     ierr = MatDenseGetLDA(B,&lda);CHKERRQ(ierr);
33*75ab9b9fSStefano Zampini     ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
34*75ab9b9fSStefano Zampini     for (j=0;j<n;j++) {
35*75ab9b9fSStefano Zampini       for (i=m;i<lda;i++) {
36*75ab9b9fSStefano Zampini         if (Bb[j*lda +i] != MAGIC_NUMBER) wBv = PETSC_TRUE;
37*75ab9b9fSStefano Zampini       }
38*75ab9b9fSStefano Zampini     }
396280154eSStefano Zampini     ierr = MatDenseRestoreArrayRead(B,&Bb);CHKERRQ(ierr);
406280154eSStefano Zampini   }
41456288a8SStefano Zampini   if (wA || wB) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong array in first Mat? %d, Wrong array in second Mat? %d",wA,wB);
42*75ab9b9fSStefano Zampini   if (wAv || wBv) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong data in first Mat? %d, Wrong data in second Mat? %d",wAv,wBv);
43*75ab9b9fSStefano Zampini   PetscFunctionReturn(0);
44*75ab9b9fSStefano Zampini }
45*75ab9b9fSStefano Zampini 
46*75ab9b9fSStefano Zampini typedef struct {
47*75ab9b9fSStefano Zampini   Mat A;
48*75ab9b9fSStefano Zampini   Mat P;
49*75ab9b9fSStefano Zampini   Mat R;
50*75ab9b9fSStefano Zampini } proj_data;
51*75ab9b9fSStefano Zampini 
52*75ab9b9fSStefano Zampini PetscErrorCode proj_destroy(void *ctx)
53*75ab9b9fSStefano Zampini {
54*75ab9b9fSStefano Zampini   proj_data      *userdata = (proj_data*)ctx;
55*75ab9b9fSStefano Zampini   PetscErrorCode ierr;
56*75ab9b9fSStefano Zampini 
57*75ab9b9fSStefano Zampini   PetscFunctionBegin;
58*75ab9b9fSStefano Zampini   if (!userdata) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing userdata");
59*75ab9b9fSStefano Zampini   ierr = MatDestroy(&userdata->A);CHKERRQ(ierr);
60*75ab9b9fSStefano Zampini   ierr = MatDestroy(&userdata->P);CHKERRQ(ierr);
61*75ab9b9fSStefano Zampini   ierr = MatDestroy(&userdata->R);CHKERRQ(ierr);
62*75ab9b9fSStefano Zampini   ierr = PetscFree(userdata);CHKERRQ(ierr);
63*75ab9b9fSStefano Zampini   PetscFunctionReturn(0);
64*75ab9b9fSStefano Zampini }
65*75ab9b9fSStefano Zampini 
66*75ab9b9fSStefano Zampini PetscErrorCode proj_mult(Mat S, Vec X, Vec Y)
67*75ab9b9fSStefano Zampini {
68*75ab9b9fSStefano Zampini   Mat            A,R,P;
69*75ab9b9fSStefano Zampini   Vec            Ax,Ay;
70*75ab9b9fSStefano Zampini   Vec            Px,Py;
71*75ab9b9fSStefano Zampini   proj_data      *userdata;
72*75ab9b9fSStefano Zampini   PetscErrorCode ierr;
73*75ab9b9fSStefano Zampini 
74*75ab9b9fSStefano Zampini   PetscFunctionBegin;
75*75ab9b9fSStefano Zampini   ierr = MatShellGetContext(S,(void**)&userdata);CHKERRQ(ierr);
76*75ab9b9fSStefano Zampini   if (!userdata) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_PLIB,"Missing userdata");
77*75ab9b9fSStefano Zampini   A = userdata->A;
78*75ab9b9fSStefano Zampini   R = userdata->R;
79*75ab9b9fSStefano Zampini   P = userdata->P;
80*75ab9b9fSStefano Zampini   if (!A) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_PLIB,"Missing matrix");
81*75ab9b9fSStefano Zampini   if (!R && !P) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_PLIB,"Missing projectors");
82*75ab9b9fSStefano Zampini   if (R && P) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_PLIB,"Both projectors");
83*75ab9b9fSStefano Zampini   ierr = MatCreateVecs(A,&Ax,&Ay);CHKERRQ(ierr);
84*75ab9b9fSStefano Zampini   if (R) {
85*75ab9b9fSStefano Zampini     ierr = MatCreateVecs(R,&Py,&Px);CHKERRQ(ierr);
86*75ab9b9fSStefano Zampini   } else {
87*75ab9b9fSStefano Zampini     ierr = MatCreateVecs(P,&Px,&Py);CHKERRQ(ierr);
88*75ab9b9fSStefano Zampini   }
89*75ab9b9fSStefano Zampini   ierr = VecCopy(X,Px);CHKERRQ(ierr);
90*75ab9b9fSStefano Zampini   if (P) {
91*75ab9b9fSStefano Zampini     ierr = MatMult(P,Px,Py);CHKERRQ(ierr);
92*75ab9b9fSStefano Zampini   } else {
93*75ab9b9fSStefano Zampini     ierr = MatMultTranspose(R,Px,Py);CHKERRQ(ierr);
94*75ab9b9fSStefano Zampini   }
95*75ab9b9fSStefano Zampini   ierr = VecCopy(Py,Ax);CHKERRQ(ierr);
96*75ab9b9fSStefano Zampini   ierr = MatMult(A,Ax,Ay);CHKERRQ(ierr);
97*75ab9b9fSStefano Zampini   ierr = VecCopy(Ay,Py);CHKERRQ(ierr);
98*75ab9b9fSStefano Zampini   if (P) {
99*75ab9b9fSStefano Zampini     ierr = MatMultTranspose(P,Py,Px);CHKERRQ(ierr);
100*75ab9b9fSStefano Zampini   } else {
101*75ab9b9fSStefano Zampini     ierr = MatMult(R,Py,Px);CHKERRQ(ierr);
102*75ab9b9fSStefano Zampini   }
103*75ab9b9fSStefano Zampini   ierr = VecCopy(Px,Y);CHKERRQ(ierr);
104*75ab9b9fSStefano Zampini   ierr = VecDestroy(&Px);CHKERRQ(ierr);
105*75ab9b9fSStefano Zampini   ierr = VecDestroy(&Py);CHKERRQ(ierr);
106*75ab9b9fSStefano Zampini   ierr = VecDestroy(&Ax);CHKERRQ(ierr);
107*75ab9b9fSStefano Zampini   ierr = VecDestroy(&Ay);CHKERRQ(ierr);
108*75ab9b9fSStefano Zampini   PetscFunctionReturn(0);
109*75ab9b9fSStefano Zampini }
110*75ab9b9fSStefano Zampini 
111*75ab9b9fSStefano Zampini PetscErrorCode MyPtShellPMultSymbolic(Mat S, Mat P, Mat PtAP, void** ctx)
112*75ab9b9fSStefano Zampini {
113*75ab9b9fSStefano Zampini   PetscErrorCode ierr;
114*75ab9b9fSStefano Zampini   proj_data      *userdata;
115*75ab9b9fSStefano Zampini 
116*75ab9b9fSStefano Zampini   PetscFunctionBegin;
117*75ab9b9fSStefano Zampini   ierr = PetscNew(&userdata);CHKERRQ(ierr);
118*75ab9b9fSStefano Zampini   ierr = MatShellSetContext(PtAP,(void*)userdata);CHKERRQ(ierr);
119*75ab9b9fSStefano Zampini   *ctx = (void *)userdata;
120*75ab9b9fSStefano Zampini   PetscFunctionReturn(0);
121*75ab9b9fSStefano Zampini }
122*75ab9b9fSStefano Zampini 
123*75ab9b9fSStefano Zampini PetscErrorCode MyPtShellPMultNumeric(Mat S, Mat P, Mat PtAP, void *ctx)
124*75ab9b9fSStefano Zampini {
125*75ab9b9fSStefano Zampini   Mat            A;
126*75ab9b9fSStefano Zampini   PetscErrorCode ierr;
127*75ab9b9fSStefano Zampini   proj_data      *userdata = (proj_data*)ctx;
128*75ab9b9fSStefano Zampini 
129*75ab9b9fSStefano Zampini   PetscFunctionBegin;
130*75ab9b9fSStefano Zampini   ierr = MatShellGetContext(S,(void**)&A);CHKERRQ(ierr);
131*75ab9b9fSStefano Zampini   ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
132*75ab9b9fSStefano Zampini   ierr = PetscObjectReference((PetscObject)P);CHKERRQ(ierr);
133*75ab9b9fSStefano Zampini   ierr = MatDestroy(&userdata->A);CHKERRQ(ierr);
134*75ab9b9fSStefano Zampini   ierr = MatDestroy(&userdata->P);CHKERRQ(ierr);
135*75ab9b9fSStefano Zampini   ierr = MatDestroy(&userdata->R);CHKERRQ(ierr);
136*75ab9b9fSStefano Zampini   userdata->A = A;
137*75ab9b9fSStefano Zampini   userdata->P = P;
138*75ab9b9fSStefano Zampini   ierr = MatShellSetOperation(PtAP,MATOP_MULT,(void (*)(void))proj_mult);CHKERRQ(ierr);
139*75ab9b9fSStefano Zampini   ierr = MatSetUp(PtAP);CHKERRQ(ierr);
140*75ab9b9fSStefano Zampini   ierr = MatAssemblyBegin(PtAP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
141*75ab9b9fSStefano Zampini   ierr = MatAssemblyEnd(PtAP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
142*75ab9b9fSStefano Zampini   PetscFunctionReturn(0);
143*75ab9b9fSStefano Zampini }
144*75ab9b9fSStefano Zampini 
145*75ab9b9fSStefano Zampini PetscErrorCode MyRShellRtMultSymbolic(Mat S, Mat R, Mat RARt, void **ctx)
146*75ab9b9fSStefano Zampini {
147*75ab9b9fSStefano Zampini   PetscErrorCode ierr;
148*75ab9b9fSStefano Zampini   proj_data      *userdata;
149*75ab9b9fSStefano Zampini 
150*75ab9b9fSStefano Zampini   PetscFunctionBegin;
151*75ab9b9fSStefano Zampini   ierr = PetscNew(&userdata);CHKERRQ(ierr);
152*75ab9b9fSStefano Zampini   ierr = MatShellSetContext(RARt,(void*)userdata);CHKERRQ(ierr);
153*75ab9b9fSStefano Zampini   *ctx = (void *)userdata;
154*75ab9b9fSStefano Zampini   PetscFunctionReturn(0);
155*75ab9b9fSStefano Zampini }
156*75ab9b9fSStefano Zampini 
157*75ab9b9fSStefano Zampini PetscErrorCode MyRShellRtMultNumeric(Mat S, Mat R, Mat RARt, void *ctx)
158*75ab9b9fSStefano Zampini {
159*75ab9b9fSStefano Zampini   Mat            A;
160*75ab9b9fSStefano Zampini   PetscErrorCode ierr;
161*75ab9b9fSStefano Zampini   proj_data      *userdata = (proj_data*)ctx;
162*75ab9b9fSStefano Zampini 
163*75ab9b9fSStefano Zampini   PetscFunctionBegin;
164*75ab9b9fSStefano Zampini   ierr = MatShellGetContext(S,(void**)&A);CHKERRQ(ierr);
165*75ab9b9fSStefano Zampini   ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
166*75ab9b9fSStefano Zampini   ierr = PetscObjectReference((PetscObject)R);CHKERRQ(ierr);
167*75ab9b9fSStefano Zampini   ierr = MatDestroy(&userdata->A);CHKERRQ(ierr);
168*75ab9b9fSStefano Zampini   ierr = MatDestroy(&userdata->P);CHKERRQ(ierr);
169*75ab9b9fSStefano Zampini   ierr = MatDestroy(&userdata->R);CHKERRQ(ierr);
170*75ab9b9fSStefano Zampini   userdata->A = A;
171*75ab9b9fSStefano Zampini   userdata->R = R;
172*75ab9b9fSStefano Zampini   ierr = MatShellSetOperation(RARt,MATOP_MULT,(void (*)(void))proj_mult);CHKERRQ(ierr);
173*75ab9b9fSStefano Zampini   ierr = MatSetUp(RARt);CHKERRQ(ierr);
174*75ab9b9fSStefano Zampini   ierr = MatAssemblyBegin(RARt,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
175*75ab9b9fSStefano Zampini   ierr = MatAssemblyEnd(RARt,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
176*75ab9b9fSStefano Zampini   PetscFunctionReturn(0);
177*75ab9b9fSStefano Zampini }
178*75ab9b9fSStefano Zampini 
179*75ab9b9fSStefano Zampini PetscErrorCode MyMatShellMatMultNumeric(Mat S, Mat B, Mat C, void *ctx)
180*75ab9b9fSStefano Zampini {
181*75ab9b9fSStefano Zampini   PetscErrorCode ierr;
182*75ab9b9fSStefano Zampini   Mat            A;
183*75ab9b9fSStefano Zampini 
184*75ab9b9fSStefano Zampini   PetscFunctionBegin;
185*75ab9b9fSStefano Zampini   ierr = MatShellGetContext(S,(void**)&A);CHKERRQ(ierr);
186*75ab9b9fSStefano Zampini   ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
187*75ab9b9fSStefano Zampini   PetscFunctionReturn(0);
188*75ab9b9fSStefano Zampini }
189*75ab9b9fSStefano Zampini 
190*75ab9b9fSStefano Zampini PetscErrorCode MyMatTransposeShellMatMultNumeric(Mat S, Mat B, Mat C, void *ctx)
191*75ab9b9fSStefano Zampini {
192*75ab9b9fSStefano Zampini   PetscErrorCode ierr;
193*75ab9b9fSStefano Zampini   Mat            A;
194*75ab9b9fSStefano Zampini 
195*75ab9b9fSStefano Zampini   PetscFunctionBegin;
196*75ab9b9fSStefano Zampini   ierr = MatShellGetContext(S,(void**)&A);CHKERRQ(ierr);
197*75ab9b9fSStefano Zampini   ierr = MatTransposeMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
198*75ab9b9fSStefano Zampini   PetscFunctionReturn(0);
199*75ab9b9fSStefano Zampini }
200*75ab9b9fSStefano Zampini 
201*75ab9b9fSStefano Zampini PetscErrorCode MyMatShellMatTransposeMultNumeric(Mat S, Mat B, Mat C, void *ctx)
202*75ab9b9fSStefano Zampini {
203*75ab9b9fSStefano Zampini   PetscErrorCode ierr;
204*75ab9b9fSStefano Zampini   Mat            A;
205*75ab9b9fSStefano Zampini 
206*75ab9b9fSStefano Zampini   PetscFunctionBegin;
207*75ab9b9fSStefano Zampini   ierr = MatShellGetContext(S,(void**)&A);CHKERRQ(ierr);
208*75ab9b9fSStefano Zampini   ierr = MatMatTransposeMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
2096280154eSStefano Zampini   PetscFunctionReturn(0);
2106280154eSStefano Zampini }
2116280154eSStefano Zampini 
2126280154eSStefano Zampini int main(int argc,char **args)
2136280154eSStefano Zampini {
214*75ab9b9fSStefano Zampini   Mat            X,B,A,Bt,T,T2,PtAP = NULL,RARt = NULL, R = NULL;
215*75ab9b9fSStefano Zampini   Vec            r,l,rs,ls;
216*75ab9b9fSStefano Zampini   PetscInt       m,n,k,M = 10,N = 10,K = 5, ldx = 3, ldb = 5, ldr = 4;
2176280154eSStefano Zampini   const char     *deft = MATAIJ;
2186280154eSStefano Zampini   char           mattype[256];
2196280154eSStefano Zampini   PetscBool      flg,symm = PETSC_FALSE,testtt = PETSC_TRUE, testnest = PETSC_TRUE, testtranspose = PETSC_TRUE, testcircular = PETSC_FALSE, local = PETSC_TRUE;
220*75ab9b9fSStefano Zampini   PetscBool      xgpu = PETSC_FALSE, bgpu = PETSC_FALSE, testshellops = PETSC_FALSE, testproj = PETSC_TRUE, testrart = PETSC_TRUE, testmatmatt = PETSC_TRUE, testmattmat = PETSC_TRUE;
221*75ab9b9fSStefano Zampini   PetscScalar    *dataX = NULL,*dataB = NULL, *dataR = NULL, *dataBt = NULL;
222*75ab9b9fSStefano Zampini   PetscScalar    *aX,*aB,*aBt;
223456288a8SStefano Zampini   PetscReal      err;
2246280154eSStefano Zampini   PetscErrorCode ierr;
2256280154eSStefano Zampini 
2266280154eSStefano Zampini   ierr = PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr;
2276280154eSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr);
2286280154eSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
2296280154eSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL);CHKERRQ(ierr);
2306280154eSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-symm",&symm,NULL);CHKERRQ(ierr);
2316280154eSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-local",&local,NULL);CHKERRQ(ierr);
2326280154eSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-ldx",&ldx,NULL);CHKERRQ(ierr);
2336280154eSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-ldb",&ldb,NULL);CHKERRQ(ierr);
234*75ab9b9fSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-ldr",&ldr,NULL);CHKERRQ(ierr);
2356280154eSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-testtranspose",&testtranspose,NULL);CHKERRQ(ierr);
2366280154eSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-testnest",&testnest,NULL);CHKERRQ(ierr);
2376280154eSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-testtt",&testtt,NULL);CHKERRQ(ierr);
2386280154eSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-testcircular",&testcircular,NULL);CHKERRQ(ierr);
239*75ab9b9fSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-testshellops",&testshellops,NULL);CHKERRQ(ierr);
240*75ab9b9fSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-testproj",&testproj,NULL);CHKERRQ(ierr);
241*75ab9b9fSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-testrart",&testrart,NULL);CHKERRQ(ierr);
242*75ab9b9fSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-testmatmatt",&testmatmatt,NULL);CHKERRQ(ierr);
243*75ab9b9fSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-testmattmat",&testmattmat,NULL);CHKERRQ(ierr);
244456288a8SStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-xgpu",&xgpu,NULL);CHKERRQ(ierr);
245456288a8SStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-bgpu",&bgpu,NULL);CHKERRQ(ierr);
246*75ab9b9fSStefano Zampini   ierr = PetscOptionsGetScalar(NULL,NULL,"-magic_number",&MAGIC_NUMBER,NULL);CHKERRQ(ierr);
247*75ab9b9fSStefano Zampini   if (M != N) testproj = PETSC_FALSE;
248*75ab9b9fSStefano Zampini 
2496280154eSStefano Zampini   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
2506280154eSStefano Zampini   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
2516280154eSStefano Zampini   ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
2526280154eSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
2536280154eSStefano Zampini   ierr = MatSetRandom(A,NULL);CHKERRQ(ierr);
2546280154eSStefano Zampini   if (M==N && symm) {
2556280154eSStefano Zampini     Mat AT;
2566280154eSStefano Zampini 
257*75ab9b9fSStefano Zampini     ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr);
2586280154eSStefano Zampini     ierr = MatAXPY(A,1.0,AT,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2596280154eSStefano Zampini     ierr = MatDestroy(&AT);CHKERRQ(ierr);
260*75ab9b9fSStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
2616280154eSStefano Zampini   }
2626280154eSStefano Zampini   ierr = MatViewFromOptions(A,NULL,"-A_init_view");CHKERRQ(ierr);
2636280154eSStefano Zampini   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","","");CHKERRQ(ierr);
264456288a8SStefano Zampini   ierr = PetscOptionsFList("-A_mat_type","Matrix type","MatSetType",MatList,deft,mattype,256,&flg);CHKERRQ(ierr);
2656280154eSStefano Zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2666280154eSStefano Zampini   if (flg) {
2676280154eSStefano Zampini     Mat A2;
2686280154eSStefano Zampini 
269456288a8SStefano Zampini     /* MATSEQAIJCUSPARSE does not support MAT_INITIAL_MATRIX */
270456288a8SStefano Zampini     ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
271456288a8SStefano Zampini     ierr = MatConvert(A,mattype,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
2726280154eSStefano Zampini     ierr = MatMultEqual(A,A2,10,&flg);CHKERRQ(ierr);
2736280154eSStefano Zampini     if (!flg) {
2746280154eSStefano Zampini       Mat AE,A2E;
2756280154eSStefano Zampini 
2766280154eSStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with convert\n");CHKERRQ(ierr);
2776280154eSStefano Zampini       ierr = MatComputeOperator(A,MATDENSE,&AE);CHKERRQ(ierr);
2786280154eSStefano Zampini       ierr = MatComputeOperator(A2,MATDENSE,&A2E);CHKERRQ(ierr);
2796280154eSStefano Zampini       ierr = MatView(AE,NULL);CHKERRQ(ierr);
2806280154eSStefano Zampini       ierr = MatView(A2E,NULL);CHKERRQ(ierr);
2816280154eSStefano Zampini       ierr = MatAXPY(A2E,-1.0,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
2826280154eSStefano Zampini       ierr = MatView(A2E,NULL);CHKERRQ(ierr);
2836280154eSStefano Zampini       ierr = MatDestroy(&A2E);CHKERRQ(ierr);
2846280154eSStefano Zampini       ierr = MatDestroy(&AE);CHKERRQ(ierr);
2856280154eSStefano Zampini     }
286456288a8SStefano Zampini     ierr = MatDestroy(&A2);CHKERRQ(ierr);
2876280154eSStefano Zampini   }
2886280154eSStefano Zampini   ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr);
2896280154eSStefano Zampini 
2906280154eSStefano Zampini   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
2916280154eSStefano Zampini   if (local) {
292*75ab9b9fSStefano Zampini     PetscInt i;
293*75ab9b9fSStefano Zampini 
2946280154eSStefano Zampini     ierr = PetscMalloc1(PetscMax((m+ldx)*K,1),&dataX);CHKERRQ(ierr);
2956280154eSStefano Zampini     ierr = PetscMalloc1(PetscMax((n+ldb)*K,1),&dataB);CHKERRQ(ierr);
296*75ab9b9fSStefano Zampini     for (i=0;i<PetscMax((m+ldx)*K,1);i++) dataX[i] = MAGIC_NUMBER;
297*75ab9b9fSStefano Zampini     for (i=0;i<PetscMax((n+ldb)*K,1);i++) dataB[i] = MAGIC_NUMBER;
2986280154eSStefano Zampini   }
2996280154eSStefano Zampini   ierr = MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,N,K,dataB,&B);CHKERRQ(ierr);
3006280154eSStefano Zampini   ierr = MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,dataX,&X);CHKERRQ(ierr);
301*75ab9b9fSStefano Zampini   if (local) {
302*75ab9b9fSStefano Zampini     ierr = MatDenseSetLDA(X,m+ldx);CHKERRQ(ierr);
303*75ab9b9fSStefano Zampini     ierr = MatDenseSetLDA(B,n+ldb);CHKERRQ(ierr);
304*75ab9b9fSStefano Zampini   }
305*75ab9b9fSStefano Zampini   ierr = MatGetLocalSize(B,NULL,&k);CHKERRQ(ierr);
306*75ab9b9fSStefano Zampini   if (local) {
307*75ab9b9fSStefano Zampini     PetscInt i;
308*75ab9b9fSStefano Zampini 
309*75ab9b9fSStefano Zampini     ierr = PetscMalloc1(PetscMax((k+ldr)*N,1),&dataBt);CHKERRQ(ierr);
310*75ab9b9fSStefano Zampini     for (i=0;i<PetscMax((k+ldr)*N,1);i++) dataBt[i] = MAGIC_NUMBER;
311*75ab9b9fSStefano Zampini   }
312*75ab9b9fSStefano Zampini   ierr = MatCreateDense(PETSC_COMM_WORLD,k,n,K,N,dataBt,&Bt);CHKERRQ(ierr);
313*75ab9b9fSStefano Zampini   if (local) {
314*75ab9b9fSStefano Zampini     ierr = MatDenseSetLDA(Bt,k+ldr);CHKERRQ(ierr);
315*75ab9b9fSStefano Zampini   }
3166280154eSStefano Zampini 
3176280154eSStefano Zampini   /* store pointer to dense data for testing */
3186280154eSStefano Zampini   ierr = MatDenseGetArrayRead(B,(const PetscScalar**)&dataB);CHKERRQ(ierr);
3196280154eSStefano Zampini   ierr = MatDenseGetArrayRead(X,(const PetscScalar**)&dataX);CHKERRQ(ierr);
320*75ab9b9fSStefano Zampini   ierr = MatDenseGetArrayRead(Bt,(const PetscScalar**)&dataBt);CHKERRQ(ierr);
3216280154eSStefano Zampini   aX   = dataX;
3226280154eSStefano Zampini   aB   = dataB;
323*75ab9b9fSStefano Zampini   aBt  = dataBt;
324*75ab9b9fSStefano Zampini   ierr = MatDenseRestoreArrayRead(Bt,(const PetscScalar**)&dataBt);CHKERRQ(ierr);
3256280154eSStefano Zampini   ierr = MatDenseRestoreArrayRead(B,(const PetscScalar**)&dataB);CHKERRQ(ierr);
3266280154eSStefano Zampini   ierr = MatDenseRestoreArrayRead(X,(const PetscScalar**)&dataX);CHKERRQ(ierr);
3276280154eSStefano Zampini   if (local) {
3286280154eSStefano Zampini     dataX  = aX;
3296280154eSStefano Zampini     dataB  = aB;
330*75ab9b9fSStefano Zampini     dataBt = aBt;
3316280154eSStefano Zampini   }
332*75ab9b9fSStefano Zampini 
333*75ab9b9fSStefano Zampini   ierr = MatSetRandom(X,NULL);CHKERRQ(ierr);
3346280154eSStefano Zampini   ierr = MatSetRandom(B,NULL);CHKERRQ(ierr);
335*75ab9b9fSStefano Zampini   ierr = MatSetRandom(Bt,NULL);CHKERRQ(ierr);
336*75ab9b9fSStefano Zampini   ierr = CheckLocal(X,NULL,aX,NULL);CHKERRQ(ierr);
337*75ab9b9fSStefano Zampini   ierr = CheckLocal(Bt,B,aBt,aB);CHKERRQ(ierr);
338*75ab9b9fSStefano Zampini 
339456288a8SStefano Zampini   /* convert to CUDA if needed */
340456288a8SStefano Zampini   if (bgpu) {
341456288a8SStefano Zampini     ierr = MatConvert(B,MATDENSECUDA,MAT_INPLACE_MATRIX,&B);CHKERRQ(ierr);
342*75ab9b9fSStefano Zampini     ierr = MatConvert(Bt,MATDENSECUDA,MAT_INPLACE_MATRIX,&Bt);CHKERRQ(ierr);
343456288a8SStefano Zampini   }
344456288a8SStefano Zampini   if (xgpu) {
345456288a8SStefano Zampini     ierr = MatConvert(X,MATDENSECUDA,MAT_INPLACE_MATRIX,&X);CHKERRQ(ierr);
346456288a8SStefano Zampini   }
347456288a8SStefano Zampini   ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
3486280154eSStefano Zampini 
3496280154eSStefano Zampini   /* Test reusing a previously allocated dense buffer */
3506280154eSStefano Zampini   ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
3516280154eSStefano Zampini   ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
3526280154eSStefano Zampini   ierr = MatMatMultEqual(A,B,X,10,&flg);CHKERRQ(ierr);
3536280154eSStefano Zampini   if (!flg) {
3546280154eSStefano Zampini     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with reusage\n");CHKERRQ(ierr);
3556280154eSStefano Zampini     ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
3566280154eSStefano Zampini     ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
3576280154eSStefano Zampini     ierr = MatView(T,NULL);CHKERRQ(ierr);
3586280154eSStefano Zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
3596280154eSStefano Zampini   }
3606280154eSStefano Zampini 
361*75ab9b9fSStefano Zampini   /* Test MatTransposeMat and MatMatTranspose */
362*75ab9b9fSStefano Zampini   if (testmattmat) {
363*75ab9b9fSStefano Zampini     ierr = MatTransposeMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr);
364*75ab9b9fSStefano Zampini     ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
365*75ab9b9fSStefano Zampini     ierr = MatTransposeMatMultEqual(A,X,B,10,&flg);CHKERRQ(ierr);
366*75ab9b9fSStefano Zampini     if (!flg) {
367*75ab9b9fSStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MatTransposeMat)\n");CHKERRQ(ierr);
368*75ab9b9fSStefano Zampini       ierr = MatTransposeMatMult(A,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr);
369*75ab9b9fSStefano Zampini       ierr = MatAXPY(T,-1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
370*75ab9b9fSStefano Zampini       ierr = MatView(T,NULL);CHKERRQ(ierr);
371*75ab9b9fSStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
372*75ab9b9fSStefano Zampini     }
373*75ab9b9fSStefano Zampini   }
374*75ab9b9fSStefano Zampini   if (testmatmatt) {
375*75ab9b9fSStefano Zampini     ierr = MatMatTransposeMult(A,Bt,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
376*75ab9b9fSStefano Zampini     ierr = CheckLocal(Bt,X,aBt,aX);CHKERRQ(ierr);
377*75ab9b9fSStefano Zampini     ierr = MatMatTransposeMultEqual(A,Bt,X,10,&flg);CHKERRQ(ierr);
378*75ab9b9fSStefano Zampini     if (!flg) {
379*75ab9b9fSStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MatMatTranspose)\n");CHKERRQ(ierr);
380*75ab9b9fSStefano Zampini       ierr = MatMatTransposeMult(A,Bt,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
381*75ab9b9fSStefano Zampini       ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
382*75ab9b9fSStefano Zampini       ierr = MatView(T,NULL);CHKERRQ(ierr);
383*75ab9b9fSStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
384*75ab9b9fSStefano Zampini     }
385*75ab9b9fSStefano Zampini   }
386*75ab9b9fSStefano Zampini 
387*75ab9b9fSStefano Zampini   /* Test projection operations (PtAP and RARt) */
388*75ab9b9fSStefano Zampini   if (testproj) {
389*75ab9b9fSStefano Zampini     ierr = MatPtAP(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&PtAP);CHKERRQ(ierr);
390*75ab9b9fSStefano Zampini     ierr = MatPtAPMultEqual(A,B,PtAP,10,&flg);CHKERRQ(ierr);
391*75ab9b9fSStefano Zampini     if (!flg) {
392*75ab9b9fSStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with PtAP\n");CHKERRQ(ierr);
393*75ab9b9fSStefano Zampini       ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
394*75ab9b9fSStefano Zampini       ierr = MatTransposeMatMult(B,T,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T2);CHKERRQ(ierr);
395*75ab9b9fSStefano Zampini       ierr = MatAXPY(T2,-1.0,PtAP,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
396*75ab9b9fSStefano Zampini       ierr = MatView(T2,NULL);CHKERRQ(ierr);
397*75ab9b9fSStefano Zampini       ierr = MatDestroy(&T2);CHKERRQ(ierr);
398*75ab9b9fSStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
399*75ab9b9fSStefano Zampini     }
400*75ab9b9fSStefano Zampini     ierr = PetscMalloc1(PetscMax((k+ldr)*M,1),&dataR);CHKERRQ(ierr);
401*75ab9b9fSStefano Zampini     ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,m,K,M,dataR,&R);CHKERRQ(ierr);
402*75ab9b9fSStefano Zampini     ierr = MatDenseSetLDA(R,k+ldr);CHKERRQ(ierr);
403*75ab9b9fSStefano Zampini     ierr = MatSetRandom(R,NULL);CHKERRQ(ierr);
404*75ab9b9fSStefano Zampini     if (testrart) { /* fails for AIJCUSPARSE because RA operation is not defined */
405*75ab9b9fSStefano Zampini       ierr = MatRARt(A,R,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&RARt);CHKERRQ(ierr);
406*75ab9b9fSStefano Zampini       ierr = MatRARtMultEqual(A,R,RARt,10,&flg);CHKERRQ(ierr);
407*75ab9b9fSStefano Zampini       if (!flg) {
408*75ab9b9fSStefano Zampini         ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with RARt\n");CHKERRQ(ierr);
409*75ab9b9fSStefano Zampini         ierr = MatMatTransposeMult(A,R,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
410*75ab9b9fSStefano Zampini         ierr = MatMatMult(R,T,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T2);CHKERRQ(ierr);
411*75ab9b9fSStefano Zampini         ierr = MatAXPY(T2,-1.0,RARt,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
412*75ab9b9fSStefano Zampini         ierr = MatView(T2,NULL);CHKERRQ(ierr);
413*75ab9b9fSStefano Zampini         ierr = MatDestroy(&T2);CHKERRQ(ierr);
414*75ab9b9fSStefano Zampini         ierr = MatDestroy(&T);CHKERRQ(ierr);
415*75ab9b9fSStefano Zampini       }
416*75ab9b9fSStefano Zampini     }
417*75ab9b9fSStefano Zampini   }
418*75ab9b9fSStefano Zampini 
419456288a8SStefano Zampini   /* Test MatDenseGetColumnVec and friends */
420*75ab9b9fSStefano Zampini   ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
421456288a8SStefano Zampini   ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
422456288a8SStefano Zampini   ierr = MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&T2);CHKERRQ(ierr);
423456288a8SStefano Zampini   for (k=0;k<K;k++) {
424456288a8SStefano Zampini     Vec Xv,Tv,T2v;
425456288a8SStefano Zampini 
426456288a8SStefano Zampini     ierr = MatDenseGetColumnVecRead(X,k,&Xv);CHKERRQ(ierr);
427456288a8SStefano Zampini     ierr = MatDenseGetColumnVec(T,k,&Tv);CHKERRQ(ierr);
428456288a8SStefano Zampini     ierr = MatDenseGetColumnVecWrite(T2,k,&T2v);CHKERRQ(ierr);
429456288a8SStefano Zampini     ierr = VecCopy(Xv,T2v);CHKERRQ(ierr);
430456288a8SStefano Zampini     ierr = VecAXPY(Tv,-1.,Xv);CHKERRQ(ierr);
431456288a8SStefano Zampini     ierr = MatDenseRestoreColumnVecRead(X,k,&Xv);CHKERRQ(ierr);
432456288a8SStefano Zampini     ierr = MatDenseRestoreColumnVec(T,k,&Tv);CHKERRQ(ierr);
433456288a8SStefano Zampini     ierr = MatDenseRestoreColumnVecWrite(T2,k,&T2v);CHKERRQ(ierr);
434456288a8SStefano Zampini   }
435456288a8SStefano Zampini   ierr = MatNorm(T,NORM_FROBENIUS,&err);CHKERRQ(ierr);
436456288a8SStefano Zampini   if (err > PETSC_SMALL) {
437456288a8SStefano Zampini     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with MatDenseGetColumnVec\n");CHKERRQ(ierr);
438456288a8SStefano Zampini     ierr = MatView(T,NULL);CHKERRQ(ierr);
439456288a8SStefano Zampini   }
440456288a8SStefano Zampini   ierr = MatAXPY(T2,-1.,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
441456288a8SStefano Zampini   ierr = MatNorm(T2,NORM_FROBENIUS,&err);CHKERRQ(ierr);
442456288a8SStefano Zampini   if (err > PETSC_SMALL) {
443456288a8SStefano Zampini     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with MatDenseGetColumnVecWrite\n");CHKERRQ(ierr);
444456288a8SStefano Zampini     ierr = MatView(T2,NULL);CHKERRQ(ierr);
445456288a8SStefano Zampini   }
446456288a8SStefano Zampini   ierr = MatDestroy(&T);CHKERRQ(ierr);
447456288a8SStefano Zampini   ierr = MatDestroy(&T2);CHKERRQ(ierr);
448456288a8SStefano Zampini 
449456288a8SStefano Zampini   /* Test with MatShell */
450*75ab9b9fSStefano Zampini   ierr = MatDuplicate(A,MAT_COPY_VALUES,&T);CHKERRQ(ierr);
451*75ab9b9fSStefano Zampini   ierr = MatConvert(T,MATSHELL,MAT_INITIAL_MATRIX,&T2);CHKERRQ(ierr);
452*75ab9b9fSStefano Zampini   ierr = MatDestroy(&T);CHKERRQ(ierr);
453*75ab9b9fSStefano Zampini 
454*75ab9b9fSStefano Zampini   /* scale matrix */
455*75ab9b9fSStefano Zampini   ierr = MatScale(A,2.0);CHKERRQ(ierr);
456*75ab9b9fSStefano Zampini   ierr = MatScale(T2,2.0);CHKERRQ(ierr);
457*75ab9b9fSStefano Zampini   ierr = MatCreateVecs(A,&r,&l);CHKERRQ(ierr);
458*75ab9b9fSStefano Zampini   ierr = VecSetRandom(r,NULL);CHKERRQ(ierr);
459*75ab9b9fSStefano Zampini   ierr = VecSetRandom(l,NULL);CHKERRQ(ierr);
460*75ab9b9fSStefano Zampini   ierr = MatCreateVecs(T2,&rs,&ls);CHKERRQ(ierr);
461*75ab9b9fSStefano Zampini   ierr = VecCopy(r,rs);CHKERRQ(ierr);
462*75ab9b9fSStefano Zampini   ierr = VecCopy(l,ls);CHKERRQ(ierr);
463*75ab9b9fSStefano Zampini   if (testproj) {
464*75ab9b9fSStefano Zampini     ierr = MatDiagonalScale(A,r,r);CHKERRQ(ierr);
465*75ab9b9fSStefano Zampini     ierr = MatDiagonalScale(T2,rs,rs);CHKERRQ(ierr);
466*75ab9b9fSStefano Zampini   } else {
467*75ab9b9fSStefano Zampini     ierr = MatDiagonalScale(A,l,r);CHKERRQ(ierr);
468*75ab9b9fSStefano Zampini     ierr = MatDiagonalScale(T2,ls,rs);CHKERRQ(ierr);
469*75ab9b9fSStefano Zampini   }
470*75ab9b9fSStefano Zampini   ierr = MatDuplicate(A,MAT_COPY_VALUES,&T);CHKERRQ(ierr);
471*75ab9b9fSStefano Zampini   ierr = MatAXPY(A,4.5,T,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
472*75ab9b9fSStefano Zampini   ierr = MatAXPY(T2,4.5,T,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
473*75ab9b9fSStefano Zampini   ierr = MatMultEqual(T2,A,10,&flg);CHKERRQ(ierr);
474*75ab9b9fSStefano Zampini   if (!flg) {
475*75ab9b9fSStefano Zampini     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with MATSHELL (MatMult)\n");CHKERRQ(ierr);
476*75ab9b9fSStefano Zampini   }
477*75ab9b9fSStefano Zampini   ierr = MatMultTransposeEqual(T2,A,10,&flg);CHKERRQ(ierr);
478*75ab9b9fSStefano Zampini   if (!flg) {
479*75ab9b9fSStefano Zampini     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with MATSHELL (MatMultTranspose)\n");CHKERRQ(ierr);
480*75ab9b9fSStefano Zampini   }
481*75ab9b9fSStefano Zampini   ierr = MatDestroy(&T);CHKERRQ(ierr);
482*75ab9b9fSStefano Zampini   ierr = VecDestroy(&ls);CHKERRQ(ierr);
483*75ab9b9fSStefano Zampini   ierr = VecDestroy(&rs);CHKERRQ(ierr);
484*75ab9b9fSStefano Zampini   ierr = VecDestroy(&l);CHKERRQ(ierr);
485*75ab9b9fSStefano Zampini   ierr = VecDestroy(&r);CHKERRQ(ierr);
486*75ab9b9fSStefano Zampini 
487*75ab9b9fSStefano Zampini   /* recompute projections, test reusage */
488*75ab9b9fSStefano Zampini   if (PtAP) { ierr = MatPtAP(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&PtAP);CHKERRQ(ierr); }
489*75ab9b9fSStefano Zampini   if (RARt) { ierr = MatRARt(A,R,MAT_REUSE_MATRIX,PETSC_DEFAULT,&RARt);CHKERRQ(ierr); }
490*75ab9b9fSStefano Zampini   if (testshellops) { /* test callbacks for user defined MatProducts */
491*75ab9b9fSStefano Zampini     ierr = MatShellSetMatProductOperation(T2,MATPRODUCT_AB,NULL,MyMatShellMatMultNumeric,NULL,MATDENSE,MATDENSE);CHKERRQ(ierr);
492*75ab9b9fSStefano Zampini     ierr = MatShellSetMatProductOperation(T2,MATPRODUCT_AB,NULL,MyMatShellMatMultNumeric,NULL,MATDENSECUDA,MATDENSECUDA);CHKERRQ(ierr);
493*75ab9b9fSStefano Zampini     ierr = MatShellSetMatProductOperation(T2,MATPRODUCT_AtB,NULL,MyMatTransposeShellMatMultNumeric,NULL,MATDENSE,MATDENSE);CHKERRQ(ierr);
494*75ab9b9fSStefano Zampini     ierr = MatShellSetMatProductOperation(T2,MATPRODUCT_AtB,NULL,MyMatTransposeShellMatMultNumeric,NULL,MATDENSECUDA,MATDENSECUDA);CHKERRQ(ierr);
495*75ab9b9fSStefano Zampini     ierr = MatShellSetMatProductOperation(T2,MATPRODUCT_ABt,NULL,MyMatShellMatTransposeMultNumeric,NULL,MATDENSE,MATDENSE);CHKERRQ(ierr);
496*75ab9b9fSStefano Zampini     ierr = MatShellSetMatProductOperation(T2,MATPRODUCT_ABt,NULL,MyMatShellMatTransposeMultNumeric,NULL,MATDENSECUDA,MATDENSECUDA);CHKERRQ(ierr);
497*75ab9b9fSStefano Zampini     if (testproj) {
498*75ab9b9fSStefano Zampini       ierr = MatShellSetMatProductOperation(T2,MATPRODUCT_PtAP,MyPtShellPMultSymbolic,MyPtShellPMultNumeric,proj_destroy,MATDENSE,MATSHELL);CHKERRQ(ierr);
499*75ab9b9fSStefano Zampini       ierr = MatShellSetMatProductOperation(T2,MATPRODUCT_PtAP,MyPtShellPMultSymbolic,MyPtShellPMultNumeric,proj_destroy,MATDENSECUDA,MATSHELL);CHKERRQ(ierr);
500*75ab9b9fSStefano Zampini       ierr = MatShellSetMatProductOperation(T2,MATPRODUCT_RARt,MyRShellRtMultSymbolic,MyRShellRtMultNumeric,proj_destroy,MATDENSE,MATSHELL);CHKERRQ(ierr);
501*75ab9b9fSStefano Zampini       ierr = MatShellSetMatProductOperation(T2,MATPRODUCT_RARt,MyRShellRtMultSymbolic,MyRShellRtMultNumeric,proj_destroy,MATDENSECUDA,MATSHELL);CHKERRQ(ierr);
502*75ab9b9fSStefano Zampini     }
503*75ab9b9fSStefano Zampini   }
504*75ab9b9fSStefano Zampini   ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
505*75ab9b9fSStefano Zampini   /* we either use the shell operations or the loop over columns code, applying the operator */
506456288a8SStefano Zampini   ierr = MatMatMult(T2,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
507456288a8SStefano Zampini   ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
508456288a8SStefano Zampini   ierr = MatMatMultEqual(T2,B,X,10,&flg);CHKERRQ(ierr);
509456288a8SStefano Zampini   if (!flg) {
510456288a8SStefano Zampini     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MATSHELL)\n");CHKERRQ(ierr);
511456288a8SStefano Zampini     ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
512456288a8SStefano Zampini     ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
513456288a8SStefano Zampini     ierr = MatView(T,NULL);CHKERRQ(ierr);
514456288a8SStefano Zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
515456288a8SStefano Zampini   }
516*75ab9b9fSStefano Zampini   if (testproj) {
517*75ab9b9fSStefano Zampini     ierr = MatPtAPMultEqual(T2,B,PtAP,10,&flg);CHKERRQ(ierr);
518*75ab9b9fSStefano Zampini     if (!flg) {
519*75ab9b9fSStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with PtAP\n");CHKERRQ(ierr);
520*75ab9b9fSStefano Zampini     }
521*75ab9b9fSStefano Zampini     if (testshellops) { /* projections fail if the product operations are not specified */
522*75ab9b9fSStefano Zampini       ierr = MatPtAP(T2,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
523*75ab9b9fSStefano Zampini       ierr = MatPtAP(T2,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
524*75ab9b9fSStefano Zampini       ierr = MatPtAPMultEqual(T2,B,T,10,&flg);CHKERRQ(ierr);
525*75ab9b9fSStefano Zampini       if (!flg) {
526*75ab9b9fSStefano Zampini         Mat TE;
5276718818eSStefano Zampini 
528*75ab9b9fSStefano Zampini         ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with PtAP (user defined)\n");CHKERRQ(ierr);
529*75ab9b9fSStefano Zampini         ierr = MatComputeOperator(T,MATDENSE,&TE);CHKERRQ(ierr);
530*75ab9b9fSStefano Zampini         ierr = MatView(TE,NULL);CHKERRQ(ierr);
531*75ab9b9fSStefano Zampini         ierr = MatView(PtAP,NULL);CHKERRQ(ierr);
532*75ab9b9fSStefano Zampini         ierr = MatAXPY(TE,-1.0,PtAP,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
533*75ab9b9fSStefano Zampini         ierr = MatView(TE,NULL);CHKERRQ(ierr);
534*75ab9b9fSStefano Zampini         ierr = MatDestroy(&TE);CHKERRQ(ierr);
535*75ab9b9fSStefano Zampini       }
536*75ab9b9fSStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
537*75ab9b9fSStefano Zampini     }
538*75ab9b9fSStefano Zampini     if (RARt) {
539*75ab9b9fSStefano Zampini       ierr = MatRARtMultEqual(T2,R,RARt,10,&flg);CHKERRQ(ierr);
540*75ab9b9fSStefano Zampini       if (!flg) {
541*75ab9b9fSStefano Zampini         ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with RARt\n");CHKERRQ(ierr);
542*75ab9b9fSStefano Zampini       }
543*75ab9b9fSStefano Zampini     }
544*75ab9b9fSStefano Zampini     if (testshellops) {
545*75ab9b9fSStefano Zampini       ierr = MatRARt(T2,R,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
546*75ab9b9fSStefano Zampini       ierr = MatRARt(T2,R,MAT_REUSE_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
547*75ab9b9fSStefano Zampini       ierr = MatRARtMultEqual(T2,R,T,10,&flg);CHKERRQ(ierr);
548*75ab9b9fSStefano Zampini       if (!flg) {
549*75ab9b9fSStefano Zampini         Mat TE;
550*75ab9b9fSStefano Zampini 
551*75ab9b9fSStefano Zampini         ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with RARt (user defined)\n");CHKERRQ(ierr);
552*75ab9b9fSStefano Zampini         ierr = MatComputeOperator(T,MATDENSE,&TE);CHKERRQ(ierr);
553*75ab9b9fSStefano Zampini         ierr = MatView(TE,NULL);CHKERRQ(ierr);
554*75ab9b9fSStefano Zampini         if (RARt) {
555*75ab9b9fSStefano Zampini           ierr = MatView(RARt,NULL);CHKERRQ(ierr);
556*75ab9b9fSStefano Zampini           ierr = MatAXPY(TE,-1.0,RARt,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
557*75ab9b9fSStefano Zampini           ierr = MatView(TE,NULL);CHKERRQ(ierr);
558*75ab9b9fSStefano Zampini         }
559*75ab9b9fSStefano Zampini         ierr = MatDestroy(&TE);CHKERRQ(ierr);
560*75ab9b9fSStefano Zampini       }
561*75ab9b9fSStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
562*75ab9b9fSStefano Zampini     }
563*75ab9b9fSStefano Zampini   }
564*75ab9b9fSStefano Zampini 
565*75ab9b9fSStefano Zampini   if (testmattmat) { /* we either use the shell operations or the loop over columns code applying the transposed operator */
566456288a8SStefano Zampini     ierr = MatTransposeMatMult(T2,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr);
567456288a8SStefano Zampini     ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
568456288a8SStefano Zampini     ierr = MatTransposeMatMultEqual(T2,X,B,10,&flg);CHKERRQ(ierr);
569456288a8SStefano Zampini     if (!flg) {
570456288a8SStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MatTranspose, MATSHELL)\n");CHKERRQ(ierr);
5716718818eSStefano Zampini       ierr = MatTransposeMatMult(A,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
5726718818eSStefano Zampini       ierr = MatAXPY(T,-1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
5736718818eSStefano Zampini       ierr = MatView(T,NULL);CHKERRQ(ierr);
5746718818eSStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
5756718818eSStefano Zampini     }
576456288a8SStefano Zampini   }
577*75ab9b9fSStefano Zampini   if (testmatmatt && testshellops) { /* only when shell operations are set */
578*75ab9b9fSStefano Zampini     ierr = MatMatTransposeMult(T2,Bt,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
579*75ab9b9fSStefano Zampini     ierr = CheckLocal(Bt,X,aBt,aX);CHKERRQ(ierr);
580*75ab9b9fSStefano Zampini     ierr = MatMatTransposeMultEqual(T2,Bt,X,10,&flg);CHKERRQ(ierr);
581*75ab9b9fSStefano Zampini     if (!flg) {
582*75ab9b9fSStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MatMatTranspose, MATSHELL)\n");CHKERRQ(ierr);
583*75ab9b9fSStefano Zampini       ierr = MatMatTransposeMult(A,Bt,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
584*75ab9b9fSStefano Zampini       ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
585*75ab9b9fSStefano Zampini       ierr = MatView(T,NULL);CHKERRQ(ierr);
586*75ab9b9fSStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
587*75ab9b9fSStefano Zampini     }
588*75ab9b9fSStefano Zampini   }
589456288a8SStefano Zampini   ierr = MatDestroy(&T2);CHKERRQ(ierr);
590456288a8SStefano Zampini 
5916280154eSStefano Zampini   if (testnest) { /* test with MatNest */
5926280154eSStefano Zampini     Mat        NA;
593456288a8SStefano Zampini     const char *vtype;
5946280154eSStefano Zampini 
5956280154eSStefano Zampini     ierr = MatCreateNest(PETSC_COMM_WORLD,1,NULL,1,NULL,&A,&NA);CHKERRQ(ierr);
596456288a8SStefano Zampini     /* needed to test against CUSPARSE matrices */
597456288a8SStefano Zampini     ierr = MatGetVecType(A,&vtype);CHKERRQ(ierr);
598456288a8SStefano Zampini     ierr = MatSetVecType(NA,vtype);CHKERRQ(ierr);
5996280154eSStefano Zampini     ierr = MatViewFromOptions(NA,NULL,"-NA_view");CHKERRQ(ierr);
6006280154eSStefano Zampini     ierr = MatMatMult(NA,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
6016280154eSStefano Zampini     ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
6026280154eSStefano Zampini     ierr = MatMatMultEqual(NA,B,X,10,&flg);CHKERRQ(ierr);
6036280154eSStefano Zampini     if (!flg) {
6046280154eSStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with Nest\n");CHKERRQ(ierr);
6056280154eSStefano Zampini       ierr = MatMatMult(NA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
6066280154eSStefano Zampini       ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
6076280154eSStefano Zampini       ierr = MatView(T,NULL);CHKERRQ(ierr);
6086280154eSStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
6096280154eSStefano Zampini     }
6106280154eSStefano Zampini     ierr = MatDestroy(&NA);CHKERRQ(ierr);
6116280154eSStefano Zampini   }
6126280154eSStefano Zampini 
6136280154eSStefano Zampini   if (testtranspose) { /* test with Transpose */
6146280154eSStefano Zampini     Mat TA;
6156280154eSStefano Zampini 
616*75ab9b9fSStefano Zampini     ierr = MatCreateTranspose(A,&TA);CHKERRQ(ierr);
6176280154eSStefano Zampini     ierr = MatMatMult(TA,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr);
6186280154eSStefano Zampini     ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
6196280154eSStefano Zampini     ierr = MatMatMultEqual(TA,X,B,10,&flg);CHKERRQ(ierr);
6206280154eSStefano Zampini     if (!flg) {
6216280154eSStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with Transpose\n");CHKERRQ(ierr);
622*75ab9b9fSStefano Zampini       ierr = MatMatMult(TA,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
623*75ab9b9fSStefano Zampini       ierr = MatAXPY(T,-1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
6246280154eSStefano Zampini       ierr = MatView(T,NULL);CHKERRQ(ierr);
6256280154eSStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
6266280154eSStefano Zampini     }
6276280154eSStefano Zampini     ierr = MatDestroy(&TA);CHKERRQ(ierr);
6286280154eSStefano Zampini   }
6296280154eSStefano Zampini 
6306280154eSStefano Zampini   if (testtt) { /* test with Transpose(Transpose) */
6316280154eSStefano Zampini     Mat TA, TTA;
6326280154eSStefano Zampini 
633*75ab9b9fSStefano Zampini     ierr = MatCreateTranspose(A,&TA);CHKERRQ(ierr);
634*75ab9b9fSStefano Zampini     ierr = MatCreateTranspose(TA,&TTA);CHKERRQ(ierr);
6356280154eSStefano Zampini     ierr = MatMatMult(TTA,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
6366280154eSStefano Zampini     ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
6376280154eSStefano Zampini     ierr = MatMatMultEqual(TTA,B,X,10,&flg);CHKERRQ(ierr);
6386280154eSStefano Zampini     if (!flg) {
6396280154eSStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with Transpose(Transpose)\n");CHKERRQ(ierr);
6406280154eSStefano Zampini       ierr = MatMatMult(TTA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
6416280154eSStefano Zampini       ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
6426280154eSStefano Zampini       ierr = MatView(T,NULL);CHKERRQ(ierr);
6436280154eSStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
6446280154eSStefano Zampini     }
6456280154eSStefano Zampini     ierr = MatDestroy(&TA);CHKERRQ(ierr);
6466280154eSStefano Zampini     ierr = MatDestroy(&TTA);CHKERRQ(ierr);
6476280154eSStefano Zampini   }
6486280154eSStefano Zampini 
6496280154eSStefano Zampini   if (testcircular) { /* test circular */
6506280154eSStefano Zampini     Mat AB;
6516280154eSStefano Zampini 
6526280154eSStefano Zampini     ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AB);CHKERRQ(ierr);
6536280154eSStefano Zampini     ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
6546280154eSStefano Zampini     ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
6556280154eSStefano Zampini     if (M == N && N == K) {
6566280154eSStefano Zampini       ierr = MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr);
6576280154eSStefano Zampini     } else {
6586280154eSStefano Zampini       ierr = MatTransposeMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr);
6596280154eSStefano Zampini     }
6606280154eSStefano Zampini     ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr);
6616280154eSStefano Zampini     ierr = MatDestroy(&AB);CHKERRQ(ierr);
6626280154eSStefano Zampini   }
6636280154eSStefano Zampini   ierr = MatDestroy(&X);CHKERRQ(ierr);
664*75ab9b9fSStefano Zampini   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
6656280154eSStefano Zampini   ierr = MatDestroy(&B);CHKERRQ(ierr);
6666280154eSStefano Zampini   ierr = MatDestroy(&A);CHKERRQ(ierr);
667*75ab9b9fSStefano Zampini   ierr = MatDestroy(&R);CHKERRQ(ierr);
668*75ab9b9fSStefano Zampini   ierr = MatDestroy(&PtAP);CHKERRQ(ierr);
669*75ab9b9fSStefano Zampini   ierr = MatDestroy(&RARt);CHKERRQ(ierr);
6706280154eSStefano Zampini   ierr = PetscFree(dataX);CHKERRQ(ierr);
6716280154eSStefano Zampini   ierr = PetscFree(dataB);CHKERRQ(ierr);
672*75ab9b9fSStefano Zampini   ierr = PetscFree(dataR);CHKERRQ(ierr);
673*75ab9b9fSStefano Zampini   ierr = PetscFree(dataBt);CHKERRQ(ierr);
6746280154eSStefano Zampini   ierr = PetscFinalize();
6756280154eSStefano Zampini   return ierr;
6766280154eSStefano Zampini }
6776280154eSStefano Zampini 
6786280154eSStefano Zampini /*TEST
6796280154eSStefano Zampini 
6806280154eSStefano Zampini   test:
6816280154eSStefano Zampini     suffix: 1
682*75ab9b9fSStefano Zampini     args: -local {{0 1}} -testshellops
6836280154eSStefano Zampini 
6846280154eSStefano Zampini   test:
6856280154eSStefano Zampini     output_file: output/ex70_1.out
686456288a8SStefano Zampini     requires: cuda
687456288a8SStefano Zampini     suffix: 1_cuda
688*75ab9b9fSStefano Zampini     args: -local {{0 1}} -xgpu {{0 1}} -bgpu {{0 1}} -A_mat_type {{seqaijcusparse seqaij}} -testnest 0 -testshellops {{0 1}}
689456288a8SStefano Zampini 
690456288a8SStefano Zampini   test:
691456288a8SStefano Zampini     TODO: VecGetSubVector seems broken with CUDA
692456288a8SStefano Zampini     output_file: output/ex70_1.out
693456288a8SStefano Zampini     requires: cuda
694456288a8SStefano Zampini     suffix: 1_cuda_broken
695456288a8SStefano Zampini     args: -local {{0 1}} -xgpu {{0 1}} -bgpu {{0 1}} -A_mat_type seqaijcusparse -testnest
696456288a8SStefano Zampini 
697456288a8SStefano Zampini   test:
698456288a8SStefano Zampini     output_file: output/ex70_1.out
6996280154eSStefano Zampini     nsize: 2
7006280154eSStefano Zampini     suffix: 1_par
701*75ab9b9fSStefano Zampini     args: -local {{0 1}} -testmatmatt 0
7026280154eSStefano Zampini 
7036280154eSStefano Zampini   test:
704456288a8SStefano Zampini     output_file: output/ex70_1.out
705456288a8SStefano Zampini     requires: cuda
706456288a8SStefano Zampini     nsize: 2
707456288a8SStefano Zampini     suffix: 1_par_cuda
708*75ab9b9fSStefano Zampini     args: -local {{0 1}} -xgpu {{0 1}} -bgpu {{0 1}} -A_mat_type {{mpiaijcusparse mpiaij}} -testnest 0 -testmatmatt 0 -matmatmult_Bbn 3
709456288a8SStefano Zampini 
710456288a8SStefano Zampini   test:
7116280154eSStefano Zampini     output_file: output/ex70_1.out
7126280154eSStefano Zampini     suffix: 2
7136280154eSStefano Zampini     nsize: 1
7146280154eSStefano Zampini     args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}}
7156280154eSStefano Zampini 
7166280154eSStefano Zampini   test:
717456288a8SStefano Zampini     requires: cuda
718456288a8SStefano Zampini     output_file: output/ex70_1.out
719456288a8SStefano Zampini     suffix: 2_cuda
720456288a8SStefano Zampini     nsize: 1
721456288a8SStefano Zampini     args: -M 7 -N 9 -K 2 -local {{0 1}} -testnest 0 -A_mat_type {{seqdensecuda seqdense}} -xgpu {{0 1}} -bgpu {{0 1}}
722456288a8SStefano Zampini 
723456288a8SStefano Zampini   test:
7246280154eSStefano Zampini     output_file: output/ex70_1.out
7256280154eSStefano Zampini     suffix: 2_par
7266280154eSStefano Zampini     nsize: 2
727*75ab9b9fSStefano Zampini     args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}} -testcircular -testmatmatt 0
7286280154eSStefano Zampini 
7296280154eSStefano Zampini   test:
730456288a8SStefano Zampini     requires: cuda
731456288a8SStefano Zampini     output_file: output/ex70_1.out
732456288a8SStefano Zampini     suffix: 2_par_cuda
733456288a8SStefano Zampini     nsize: 2
734*75ab9b9fSStefano Zampini     args: -M 11 -N 9 -K 1 -local {{0 1}} -testcircular 0 -A_mat_type mpiaijcusparse -xgpu -bgpu -testnest 0 -testmatmatt 0
735456288a8SStefano Zampini 
736456288a8SStefano Zampini   test:
7376280154eSStefano Zampini     output_file: output/ex70_1.out
7386280154eSStefano Zampini     suffix: 3
739456288a8SStefano Zampini     nsize: {{1 3}}
740*75ab9b9fSStefano Zampini     args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -A_mat_type sbaij -symm -testproj 0 -testmatmatt 0
7416280154eSStefano Zampini 
7426280154eSStefano Zampini   test:
7436280154eSStefano Zampini     output_file: output/ex70_1.out
7446280154eSStefano Zampini     suffix: 4
7456280154eSStefano Zampini     nsize: 1
7466280154eSStefano Zampini     args: -M 3 -N 3 -K 3 -local {{0 1}} -testcircular
7476280154eSStefano Zampini 
7486280154eSStefano Zampini   test:
7496280154eSStefano Zampini     output_file: output/ex70_1.out
7506280154eSStefano Zampini     suffix: 5
7516280154eSStefano Zampini     nsize: {{2 4}}
752*75ab9b9fSStefano Zampini     args: -M 3 -N 3 -K 3 -local 1 -testcircular -testmatmatt 0
7536280154eSStefano Zampini 
7546280154eSStefano Zampini   test:
7556280154eSStefano Zampini     TODO: MatCreateDense broken with some processors not having local rows
7566280154eSStefano Zampini     output_file: output/ex70_1.out
7576280154eSStefano Zampini     suffix: 5_broken
7586280154eSStefano Zampini     nsize: {{2 4}}
7596280154eSStefano Zampini     args: -M 3 -N 3 -K 3 -local 0 -testcircular -testtranspose 0
7606280154eSStefano Zampini 
7616280154eSStefano Zampini   test:
7626280154eSStefano Zampini     output_file: output/ex70_1.out
7636280154eSStefano Zampini     suffix: 6
7646280154eSStefano Zampini     nsize: 1
765*75ab9b9fSStefano Zampini     args: -M {{1 3}} -N {{2 5}} -K {{1 2}} -local {{0 1}} -testcircular
7666280154eSStefano Zampini 
7676280154eSStefano Zampini   test:
7686280154eSStefano Zampini     output_file: output/ex70_1.out
7696280154eSStefano Zampini     suffix: 7
7706280154eSStefano Zampini     nsize: 1
771456288a8SStefano Zampini     args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -A_mat_type dense -testnest -testcircular
7726718818eSStefano Zampini 
7736280154eSStefano Zampini TEST*/
774