xref: /petsc/src/mat/tests/ex1.c (revision 4905a7bc61a644ac28a555b575668251734ce1fa)
1c4762a1bSJed Brown 
2*4905a7bcSToby Isaac static char help[] = "Tests LU, Cholesky, and QR factorization and MatMatSolve() for a sequential dense matrix. \n\
3c4762a1bSJed Brown                       For MATSEQDENSE matrix, the factorization is just a thin wrapper to LAPACK.       \n\
4c4762a1bSJed Brown                       For MATSEQDENSECUDA, it uses cusolverDn routines \n\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #include <petscmat.h>
7c4762a1bSJed Brown 
8*4905a7bcSToby Isaac static PetscErrorCode createMatsAndVecs(PetscInt m, PetscInt n, PetscInt nrhs, PetscBool full, Mat *_mat, Mat *_RHS, Mat *_SOLU, Vec *_x, Vec *_y, Vec *_b)
9c4762a1bSJed Brown {
10c4762a1bSJed Brown   PetscRandom    rand;
11*4905a7bcSToby Isaac   Mat            mat,RHS,SOLU;
12*4905a7bcSToby Isaac   PetscInt       rstart, rend;
13*4905a7bcSToby Isaac   PetscInt       cstart, cend;
14*4905a7bcSToby Isaac   PetscScalar    value = 1.0;
15*4905a7bcSToby Isaac   Vec            x, y, b;
16*4905a7bcSToby Isaac   PetscErrorCode ierr;
17c4762a1bSJed Brown 
18*4905a7bcSToby Isaac   PetscFunctionBegin;
19c4762a1bSJed Brown   /* create multiple vectors RHS and SOLU */
20c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&RHS);CHKERRQ(ierr);
21*4905a7bcSToby Isaac   ierr = MatSetSizes(RHS,PETSC_DECIDE,PETSC_DECIDE,m,nrhs);CHKERRQ(ierr);
22c4762a1bSJed Brown   ierr = MatSetType(RHS,MATDENSE);CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = MatSetOptionsPrefix(RHS,"rhs_");CHKERRQ(ierr);
24c4762a1bSJed Brown   ierr = MatSetFromOptions(RHS);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = MatSeqDenseSetPreallocation(RHS,NULL);CHKERRQ(ierr);
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = MatSetRandom(RHS,rand);CHKERRQ(ierr);
30c4762a1bSJed Brown 
31*4905a7bcSToby Isaac   if (m == n) {
32c4762a1bSJed Brown     ierr = MatDuplicate(RHS,MAT_DO_NOT_COPY_VALUES,&SOLU);CHKERRQ(ierr);
33*4905a7bcSToby Isaac   } else {
34*4905a7bcSToby Isaac     ierr = MatCreate(PETSC_COMM_WORLD,&SOLU);CHKERRQ(ierr);
35*4905a7bcSToby Isaac     ierr = MatSetSizes(SOLU,PETSC_DECIDE,PETSC_DECIDE,n,nrhs);CHKERRQ(ierr);
36*4905a7bcSToby Isaac     ierr = MatSetType(SOLU,MATDENSE);CHKERRQ(ierr);
37*4905a7bcSToby Isaac     ierr = MatSeqDenseSetPreallocation(SOLU,NULL);CHKERRQ(ierr);
38*4905a7bcSToby Isaac   }
39*4905a7bcSToby Isaac   ierr = MatSetRandom(SOLU,rand);CHKERRQ(ierr);
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* create matrix */
42c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&mat);CHKERRQ(ierr);
43*4905a7bcSToby Isaac   ierr = MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = MatSetType(mat,MATDENSE);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = MatSetUp(mat);CHKERRQ(ierr);
47c4762a1bSJed Brown   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
48*4905a7bcSToby Isaac   ierr = MatGetOwnershipRangeColumn(mat,&cstart,&cend);CHKERRQ(ierr);
49c4762a1bSJed Brown   if (!full) {
50*4905a7bcSToby Isaac     for (PetscInt i=rstart; i<rend; i++) {
51*4905a7bcSToby Isaac       if (m == n) {
52c4762a1bSJed Brown         value = (PetscReal)i+1;
53c4762a1bSJed Brown         ierr  = MatSetValues(mat,1,&i,1,&i,&value,INSERT_VALUES);CHKERRQ(ierr);
54*4905a7bcSToby Isaac       } else {
55*4905a7bcSToby Isaac         for (PetscInt j = cstart; j < cend; j++) {
56*4905a7bcSToby Isaac           value = ((PetscScalar)i+1.)/(PetscSqr(i - j) + 1.);
57*4905a7bcSToby Isaac           ierr  = MatSetValues(mat,1,&i,1,&j,&value,INSERT_VALUES);CHKERRQ(ierr);
58*4905a7bcSToby Isaac         }
59*4905a7bcSToby Isaac       }
60c4762a1bSJed Brown     }
61c4762a1bSJed Brown     ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
62c4762a1bSJed Brown     ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
63c4762a1bSJed Brown   } else {
64*4905a7bcSToby Isaac     ierr = MatSetRandom(mat,rand);CHKERRQ(ierr);
65*4905a7bcSToby Isaac     if (m == n) {
66c4762a1bSJed Brown       Mat T;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown       ierr = MatMatTransposeMult(mat,mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr);
69c4762a1bSJed Brown       ierr = MatDestroy(&mat);CHKERRQ(ierr);
70c4762a1bSJed Brown       mat  = T;
71c4762a1bSJed Brown     }
72*4905a7bcSToby Isaac   }
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* create single vectors */
75c4762a1bSJed Brown   ierr  = MatCreateVecs(mat,&x,&b);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr  = VecDuplicate(x,&y);CHKERRQ(ierr);
77c4762a1bSJed Brown   ierr  = VecSet(x,value);CHKERRQ(ierr);
78*4905a7bcSToby Isaac   ierr  = PetscRandomDestroy(&rand);CHKERRQ(ierr);
79*4905a7bcSToby Isaac   *_mat  = mat;
80*4905a7bcSToby Isaac   *_RHS  = RHS;
81*4905a7bcSToby Isaac   *_SOLU = SOLU;
82*4905a7bcSToby Isaac   *_x    = x;
83*4905a7bcSToby Isaac   *_y    = y;
84*4905a7bcSToby Isaac   *_b    = b;
85*4905a7bcSToby Isaac   PetscFunctionReturn(0);
86*4905a7bcSToby Isaac }
87*4905a7bcSToby Isaac 
88*4905a7bcSToby Isaac int main(int argc,char **argv)
89*4905a7bcSToby Isaac {
90*4905a7bcSToby Isaac   Mat            mat,F,RHS,SOLU;
91*4905a7bcSToby Isaac   MatInfo        info;
92*4905a7bcSToby Isaac   PetscErrorCode ierr;
93*4905a7bcSToby Isaac   PetscInt       m = 15, n = 10,i,j,nrhs=2;
94*4905a7bcSToby Isaac   Vec            x,y,b,ytmp;
95*4905a7bcSToby Isaac   IS             perm;
96*4905a7bcSToby Isaac   PetscReal      norm,tol=PETSC_SMALL;
97*4905a7bcSToby Isaac   PetscMPIInt    size;
98*4905a7bcSToby Isaac   char           solver[64];
99*4905a7bcSToby Isaac   PetscBool      inplace,full = PETSC_FALSE, ldl = PETSC_TRUE;
100*4905a7bcSToby Isaac 
101*4905a7bcSToby Isaac   ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr;
102*4905a7bcSToby Isaac   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
103*4905a7bcSToby Isaac   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
104*4905a7bcSToby Isaac   ierr = PetscStrcpy(solver,"petsc");CHKERRQ(ierr);
105*4905a7bcSToby Isaac   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
106*4905a7bcSToby Isaac   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
107*4905a7bcSToby Isaac   ierr = PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL);CHKERRQ(ierr);
108*4905a7bcSToby Isaac   ierr = PetscOptionsGetBool(NULL,NULL,"-ldl",&ldl,NULL);CHKERRQ(ierr);
109*4905a7bcSToby Isaac   ierr = PetscOptionsGetBool(NULL,NULL,"-full",&full,NULL);CHKERRQ(ierr);
110*4905a7bcSToby Isaac   ierr = PetscOptionsGetString(NULL,NULL,"-solver_type",solver,sizeof(solver),NULL);CHKERRQ(ierr);
111*4905a7bcSToby Isaac 
112*4905a7bcSToby Isaac   ierr = createMatsAndVecs(n, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b);CHKERRQ(ierr);
113*4905a7bcSToby Isaac   ierr = VecDuplicate(y,&ytmp);CHKERRQ(ierr);
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* Only SeqDense* support in-place factorizations and NULL permutations */
116c4762a1bSJed Brown   ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATSEQDENSE,&inplace);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = MatGetLocalSize(mat,&i,NULL);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = MatGetOwnershipRange(mat,&j,NULL);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,i,j,1,&perm);CHKERRQ(ierr);
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
122c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"matrix nonzeros = %D, allocated nonzeros = %D\n",
123c4762a1bSJed Brown                      (PetscInt)info.nz_used,(PetscInt)info.nz_allocated);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = MatMult(mat,x,b);CHKERRQ(ierr);
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* Cholesky factorization - perm and factinfo are ignored by LAPACK */
127c4762a1bSJed Brown   /* in-place Cholesky */
128c4762a1bSJed Brown   if (inplace) {
129c4762a1bSJed Brown     Mat RHS2;
130c4762a1bSJed Brown 
131c4762a1bSJed Brown     ierr = MatDuplicate(mat,MAT_COPY_VALUES,&F);CHKERRQ(ierr);
132c4762a1bSJed Brown     if (!ldl) { ierr = MatSetOption(F,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr); }
133c4762a1bSJed Brown     ierr = MatCholeskyFactor(F,perm,0);CHKERRQ(ierr);
134c4762a1bSJed Brown     ierr = MatSolve(F,b,y);CHKERRQ(ierr);
135c4762a1bSJed Brown     ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
136c4762a1bSJed Brown     ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
137c4762a1bSJed Brown     if (norm > tol) {
138c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for in-place Cholesky %g\n",(double)norm);CHKERRQ(ierr);
139c4762a1bSJed Brown     }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown     ierr = MatMatSolve(F,RHS,SOLU);CHKERRQ(ierr);
142c4762a1bSJed Brown     ierr = MatMatMult(mat,SOLU,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&RHS2);CHKERRQ(ierr);
143c4762a1bSJed Brown     ierr = MatAXPY(RHS,-1.0,RHS2,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
144c4762a1bSJed Brown     ierr = MatNorm(RHS,NORM_FROBENIUS,&norm);CHKERRQ(ierr);
145c4762a1bSJed Brown     if (norm > tol) {
146c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: Norm of residual for in-place Cholesky (MatMatSolve) %g\n",(double)norm);CHKERRQ(ierr);
147c4762a1bSJed Brown     }
148c4762a1bSJed Brown     ierr = MatDestroy(&F);CHKERRQ(ierr);
149c4762a1bSJed Brown     ierr = MatDestroy(&RHS2);CHKERRQ(ierr);
150c4762a1bSJed Brown   }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* out-of-place Cholesky */
153c4762a1bSJed Brown   ierr = MatGetFactor(mat,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
154c4762a1bSJed Brown   if (!ldl) { ierr = MatSetOption(F,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr); }
155c4762a1bSJed Brown   ierr = MatCholeskyFactorSymbolic(F,mat,perm,0);CHKERRQ(ierr);
156c4762a1bSJed Brown   ierr = MatCholeskyFactorNumeric(F,mat,0);CHKERRQ(ierr);
157c4762a1bSJed Brown   ierr = MatSolve(F,b,y);CHKERRQ(ierr);
158c4762a1bSJed Brown   ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
159c4762a1bSJed Brown   ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
160c4762a1bSJed Brown   if (norm > tol) {
161c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place Cholesky %g\n",(double)norm);CHKERRQ(ierr);
162c4762a1bSJed Brown   }
163c4762a1bSJed Brown   ierr = MatDestroy(&F);CHKERRQ(ierr);
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /* LU factorization - perms and factinfo are ignored by LAPACK */
166c4762a1bSJed Brown   i    = n-1;
167c4762a1bSJed Brown   ierr = MatZeroRows(mat,1,&i,-1.0,NULL,NULL);CHKERRQ(ierr);
168c4762a1bSJed Brown   ierr = MatMult(mat,x,b);CHKERRQ(ierr);
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* in-place LU */
171c4762a1bSJed Brown   if (inplace) {
172c4762a1bSJed Brown     Mat RHS2;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown     ierr = MatDuplicate(mat,MAT_COPY_VALUES,&F);CHKERRQ(ierr);
175c4762a1bSJed Brown     ierr = MatLUFactor(F,perm,perm,0);CHKERRQ(ierr);
176c4762a1bSJed Brown     ierr = MatSolve(F,b,y);CHKERRQ(ierr);
177c4762a1bSJed Brown     ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
178c4762a1bSJed Brown     ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
179c4762a1bSJed Brown     if (norm > tol) {
180c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for in-place LU %g\n",(double)norm);CHKERRQ(ierr);
181c4762a1bSJed Brown     }
182c4762a1bSJed Brown     ierr = MatMatSolve(F,RHS,SOLU);CHKERRQ(ierr);
183c4762a1bSJed Brown     ierr = MatMatMult(mat,SOLU,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&RHS2);CHKERRQ(ierr);
184c4762a1bSJed Brown     ierr = MatAXPY(RHS,-1.0,RHS2,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
185c4762a1bSJed Brown     ierr = MatNorm(RHS,NORM_FROBENIUS,&norm);CHKERRQ(ierr);
186c4762a1bSJed Brown     if (norm > tol) {
187c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: Norm of residual for in-place LU (MatMatSolve) %g\n",(double)norm);CHKERRQ(ierr);
188c4762a1bSJed Brown     }
189c4762a1bSJed Brown     ierr = MatDestroy(&F);CHKERRQ(ierr);
190c4762a1bSJed Brown     ierr = MatDestroy(&RHS2);CHKERRQ(ierr);
191c4762a1bSJed Brown   }
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   /* out-of-place LU */
194c4762a1bSJed Brown   ierr = MatGetFactor(mat,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
195c4762a1bSJed Brown   ierr = MatLUFactorSymbolic(F,mat,perm,perm,0);CHKERRQ(ierr);
196c4762a1bSJed Brown   ierr = MatLUFactorNumeric(F,mat,0);CHKERRQ(ierr);
197c4762a1bSJed Brown   ierr = MatSolve(F,b,y);CHKERRQ(ierr);
198c4762a1bSJed Brown   ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
199c4762a1bSJed Brown   ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
200c4762a1bSJed Brown   if (norm > tol) {
201c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place LU %g\n",(double)norm);CHKERRQ(ierr);
202c4762a1bSJed Brown   }
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* free space */
205c4762a1bSJed Brown   ierr = ISDestroy(&perm);CHKERRQ(ierr);
206c4762a1bSJed Brown   ierr = MatDestroy(&F);CHKERRQ(ierr);
207c4762a1bSJed Brown   ierr = MatDestroy(&mat);CHKERRQ(ierr);
208c4762a1bSJed Brown   ierr = MatDestroy(&RHS);CHKERRQ(ierr);
209c4762a1bSJed Brown   ierr = MatDestroy(&SOLU);CHKERRQ(ierr);
210*4905a7bcSToby Isaac   ierr = VecDestroy(&x);CHKERRQ(ierr);
211*4905a7bcSToby Isaac   ierr = VecDestroy(&b);CHKERRQ(ierr);
212*4905a7bcSToby Isaac   ierr = VecDestroy(&y);CHKERRQ(ierr);
213*4905a7bcSToby Isaac   ierr = VecDestroy(&ytmp);CHKERRQ(ierr);
214*4905a7bcSToby Isaac 
215*4905a7bcSToby Isaac   /* setup rectanglar */
216*4905a7bcSToby Isaac   ierr = createMatsAndVecs(m, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b);CHKERRQ(ierr);
217*4905a7bcSToby Isaac   ierr = VecDuplicate(y,&ytmp);CHKERRQ(ierr);
218*4905a7bcSToby Isaac 
219*4905a7bcSToby Isaac   /* QR factorization - perms and factinfo are ignored by LAPACK */
220*4905a7bcSToby Isaac   ierr = MatMult(mat,x,b);CHKERRQ(ierr);
221*4905a7bcSToby Isaac 
222*4905a7bcSToby Isaac   /* in-place QR */
223*4905a7bcSToby Isaac   if (inplace) {
224*4905a7bcSToby Isaac     Mat SOLU2;
225*4905a7bcSToby Isaac 
226*4905a7bcSToby Isaac     ierr = MatDuplicate(mat,MAT_COPY_VALUES,&F);CHKERRQ(ierr);
227*4905a7bcSToby Isaac     ierr = MatQRFactor(F,NULL,0);CHKERRQ(ierr);
228*4905a7bcSToby Isaac     ierr = MatSolve(F,b,y);CHKERRQ(ierr);
229*4905a7bcSToby Isaac     ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
230*4905a7bcSToby Isaac     ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
231*4905a7bcSToby Isaac     if (norm > tol) {
232*4905a7bcSToby Isaac       ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for in-place QR %g\n",(double)norm);CHKERRQ(ierr);
233*4905a7bcSToby Isaac     }
234*4905a7bcSToby Isaac     ierr = MatMatMult(mat,SOLU,MAT_REUSE_MATRIX,PETSC_DEFAULT,&RHS);CHKERRQ(ierr);
235*4905a7bcSToby Isaac     ierr = MatDuplicate(SOLU, MAT_DO_NOT_COPY_VALUES, &SOLU2);CHKERRQ(ierr);
236*4905a7bcSToby Isaac     ierr = MatMatSolve(F,RHS,SOLU2);CHKERRQ(ierr);
237*4905a7bcSToby Isaac     ierr = MatAXPY(SOLU2,-1.0,SOLU,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
238*4905a7bcSToby Isaac     ierr = MatNorm(SOLU2,NORM_FROBENIUS,&norm);CHKERRQ(ierr);
239*4905a7bcSToby Isaac     if (norm > tol) {
240*4905a7bcSToby Isaac       ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: Norm of error for in-place QR (MatMatSolve) %g\n",(double)norm);CHKERRQ(ierr);
241*4905a7bcSToby Isaac     }
242*4905a7bcSToby Isaac     ierr = MatDestroy(&F);CHKERRQ(ierr);
243*4905a7bcSToby Isaac     ierr = MatDestroy(&SOLU2);CHKERRQ(ierr);
244*4905a7bcSToby Isaac   }
245*4905a7bcSToby Isaac 
246*4905a7bcSToby Isaac   /* out-of-place QR */
247*4905a7bcSToby Isaac   ierr = MatGetFactor(mat,solver,MAT_FACTOR_QR,&F);CHKERRQ(ierr);
248*4905a7bcSToby Isaac   ierr = MatQRFactorSymbolic(F,mat,NULL,NULL);CHKERRQ(ierr);
249*4905a7bcSToby Isaac   ierr = MatQRFactorNumeric(F,mat,NULL);CHKERRQ(ierr);
250*4905a7bcSToby Isaac   ierr = MatSolve(F,b,y);CHKERRQ(ierr);
251*4905a7bcSToby Isaac   ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
252*4905a7bcSToby Isaac   ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
253*4905a7bcSToby Isaac   if (norm > tol) {
254*4905a7bcSToby Isaac     ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place QR %g\n",(double)norm);CHKERRQ(ierr);
255*4905a7bcSToby Isaac   }
256*4905a7bcSToby Isaac 
257*4905a7bcSToby Isaac   if (m == n) {
258*4905a7bcSToby Isaac     /* out-of-place MatSolveTranspose */
259*4905a7bcSToby Isaac     ierr = MatMultTranspose(mat,x,b);CHKERRQ(ierr);
260*4905a7bcSToby Isaac     ierr = MatSolveTranspose(F,b,y);CHKERRQ(ierr);
261*4905a7bcSToby Isaac     ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
262*4905a7bcSToby Isaac     ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
263*4905a7bcSToby Isaac     if (norm > tol) {
264*4905a7bcSToby Isaac       ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place QR %g\n",(double)norm);CHKERRQ(ierr);
265*4905a7bcSToby Isaac     }
266*4905a7bcSToby Isaac   }
267*4905a7bcSToby Isaac 
268*4905a7bcSToby Isaac   /* free space */
269*4905a7bcSToby Isaac   ierr = MatDestroy(&F);CHKERRQ(ierr);
270*4905a7bcSToby Isaac   ierr = MatDestroy(&mat);CHKERRQ(ierr);
271*4905a7bcSToby Isaac   ierr = MatDestroy(&RHS);CHKERRQ(ierr);
272*4905a7bcSToby Isaac   ierr = MatDestroy(&SOLU);CHKERRQ(ierr);
273c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
274c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
275c4762a1bSJed Brown   ierr = VecDestroy(&y);CHKERRQ(ierr);
276c4762a1bSJed Brown   ierr = VecDestroy(&ytmp);CHKERRQ(ierr);
277c4762a1bSJed Brown   ierr = PetscFinalize();
278c4762a1bSJed Brown   return ierr;
279c4762a1bSJed Brown }
280c4762a1bSJed Brown 
281c4762a1bSJed Brown 
282c4762a1bSJed Brown 
283c4762a1bSJed Brown /*TEST
284c4762a1bSJed Brown 
285c4762a1bSJed Brown    test:
286c4762a1bSJed Brown 
287c4762a1bSJed Brown    test:
288c4762a1bSJed Brown      requires: cuda
289c4762a1bSJed Brown      suffix: seqdensecuda
290c4762a1bSJed Brown      args: -mat_type seqdensecuda -rhs_mat_type seqdensecuda -ldl 0 -solver_type {{petsc cuda}}
291c4762a1bSJed Brown      output_file: output/ex1_1.out
292c4762a1bSJed Brown 
293c4762a1bSJed Brown    test:
294c4762a1bSJed Brown      requires: cuda
295c4762a1bSJed Brown      suffix: seqdensecuda_seqaijcusparse
296c4762a1bSJed Brown      args: -mat_type seqaijcusparse -rhs_mat_type seqdensecuda
297c4762a1bSJed Brown      output_file: output/ex1_2.out
298c4762a1bSJed Brown 
299c4762a1bSJed Brown    test:
300c4762a1bSJed Brown      requires: cuda viennacl
301c4762a1bSJed Brown      suffix: seqdensecuda_seqaijviennacl
302c4762a1bSJed Brown      args: -mat_type seqaijviennacl -rhs_mat_type seqdensecuda
303c4762a1bSJed Brown      output_file: output/ex1_2.out
304c4762a1bSJed Brown 
305*4905a7bcSToby Isaac    test:
306*4905a7bcSToby Isaac      suffix: 4
307*4905a7bcSToby Isaac      args: -m 10 -n 10
308*4905a7bcSToby Isaac      output_file: output/ex1_1.out
309*4905a7bcSToby Isaac 
310c4762a1bSJed Brown TEST*/
311