xref: /petsc/src/mat/tests/ex204.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
1c4762a1bSJed Brown static char help[] = "Test ViennaCL Matrix Conversions";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **args)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   PetscErrorCode ierr;
8c4762a1bSJed Brown   PetscMPIInt rank,size;
9c4762a1bSJed Brown 
10c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
11c4762a1bSJed Brown 
12ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
13ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   /* Construct a sequential ViennaCL matrix and vector */
16c4762a1bSJed Brown   if (!rank) {
17c4762a1bSJed Brown     Mat A_vcl;
18c4762a1bSJed Brown     Vec v_vcl,r_vcl;
19c4762a1bSJed Brown     PetscInt n = 17, m = 31,nz = 5,i,cnt,j;
20c4762a1bSJed Brown     const PetscScalar val = 1.0;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown     ierr = MatCreateSeqAIJViennaCL(PETSC_COMM_SELF,m,n,nz,NULL,&A_vcl);CHKERRQ(ierr);
23c4762a1bSJed Brown 
24c4762a1bSJed Brown     /* Add nz arbitrary entries per row in arbitrary columns */
25c4762a1bSJed Brown     for (i=0;i<m;++i) {
26c4762a1bSJed Brown       for (cnt = 0; cnt<nz; ++cnt) {
27c4762a1bSJed Brown         j = (19 * cnt + (7*i + 3)) % n;
28c4762a1bSJed Brown         ierr = MatSetValue(A_vcl,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES);CHKERRQ(ierr);
29c4762a1bSJed Brown       }
30c4762a1bSJed Brown     }
31c4762a1bSJed Brown     ierr = MatAssemblyBegin(A_vcl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
32c4762a1bSJed Brown     ierr = MatAssemblyEnd(A_vcl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
33c4762a1bSJed Brown 
34c4762a1bSJed Brown     ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl);CHKERRQ(ierr);
35c4762a1bSJed Brown     ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl);CHKERRQ(ierr);
36c4762a1bSJed Brown     ierr = VecSet(v_vcl,val);CHKERRQ(ierr);
37c4762a1bSJed Brown 
38c4762a1bSJed Brown     ierr = MatMult(A_vcl,v_vcl,r_vcl);CHKERRQ(ierr);
39c4762a1bSJed Brown 
40c4762a1bSJed Brown     ierr = VecDestroy(&v_vcl);CHKERRQ(ierr);
41c4762a1bSJed Brown     ierr = VecDestroy(&r_vcl);CHKERRQ(ierr);
42c4762a1bSJed Brown     ierr = MatDestroy(&A_vcl);CHKERRQ(ierr);
43c4762a1bSJed Brown   }
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Create a sequential AIJ matrix on rank 0 convert it to a new ViennaCL matrix, and apply it to a ViennaCL vector */
46c4762a1bSJed Brown   if (!rank) {
47c4762a1bSJed Brown     Mat               A,A_vcl;
48c4762a1bSJed Brown     Vec               v,r,v_vcl,r_vcl,d_vcl;
49c4762a1bSJed Brown     PetscInt          n=17,m=31,nz=5,i,cnt,j;
50c4762a1bSJed Brown     const PetscScalar val = 1.0;
51c4762a1bSJed Brown     PetscReal         dnorm;
52c4762a1bSJed Brown     const PetscReal   tol = 1e-5;
53c4762a1bSJed Brown 
54c4762a1bSJed Brown     ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,nz,NULL,&A);CHKERRQ(ierr);
55c4762a1bSJed Brown 
56c4762a1bSJed Brown     /* Add nz arbitrary entries per row in arbitrary columns */
57c4762a1bSJed Brown     for (i=0;i<m;++i) {
58c4762a1bSJed Brown       for (cnt = 0; cnt<nz; ++cnt) {
59c4762a1bSJed Brown         j = (19 * cnt + (7*i + 3)) % n;
60c4762a1bSJed Brown         ierr = MatSetValue(A,i,j,(PetscScalar) (0.3 * i + j),INSERT_VALUES);CHKERRQ(ierr);
61c4762a1bSJed Brown       }
62c4762a1bSJed Brown     }
63c4762a1bSJed Brown     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64c4762a1bSJed Brown     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65c4762a1bSJed Brown 
66c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix");CHKERRQ(ierr);
67c4762a1bSJed Brown 
68c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF,n,&v);CHKERRQ(ierr);
69c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF,m,&r);CHKERRQ(ierr);
70c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)r,"CPU result vector");CHKERRQ(ierr);
71c4762a1bSJed Brown     ierr = VecSet(v,val);CHKERRQ(ierr);
72c4762a1bSJed Brown     ierr = MatMult(A,v,r);CHKERRQ(ierr);
73c4762a1bSJed Brown 
74c4762a1bSJed Brown     ierr = MatConvert(A,MATSEQAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl);CHKERRQ(ierr);
75c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)A_vcl,"New ViennaCL Matrix");CHKERRQ(ierr);
76c4762a1bSJed Brown 
77c4762a1bSJed Brown     ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl);CHKERRQ(ierr);
78c4762a1bSJed Brown     ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl);CHKERRQ(ierr);
79c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector");CHKERRQ(ierr);
80c4762a1bSJed Brown     ierr = VecSet(v_vcl,val);CHKERRQ(ierr);
81c4762a1bSJed Brown     ierr = MatMult(A_vcl,v_vcl,r_vcl);CHKERRQ(ierr);
82c4762a1bSJed Brown 
83c4762a1bSJed Brown     ierr = VecDuplicate(r_vcl,&d_vcl);CHKERRQ(ierr);
84c4762a1bSJed Brown     ierr = VecCopy(r_vcl,d_vcl);CHKERRQ(ierr);
85*1e1ea65dSPierre Jolivet     ierr = VecAXPY(d_vcl,-1.0,r_vcl);CHKERRQ(ierr);
86c4762a1bSJed Brown     ierr = VecNorm(d_vcl,NORM_INFINITY,&dnorm);CHKERRQ(ierr);
87c4762a1bSJed Brown     if (dnorm > tol) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"Sequential CPU and MPI ViennaCL vector results incompatible with inf norm greater than tolerance of %g",tol);
88c4762a1bSJed Brown 
89c4762a1bSJed Brown     ierr = VecDestroy(&v);CHKERRQ(ierr);
90c4762a1bSJed Brown     ierr = VecDestroy(&r);CHKERRQ(ierr);
91c4762a1bSJed Brown     ierr = VecDestroy(&v_vcl);CHKERRQ(ierr);
92c4762a1bSJed Brown     ierr = VecDestroy(&r_vcl);CHKERRQ(ierr);
93c4762a1bSJed Brown     ierr = VecDestroy(&d_vcl);CHKERRQ(ierr);
94c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
95c4762a1bSJed Brown     ierr = MatDestroy(&A_vcl);CHKERRQ(ierr);
96c4762a1bSJed Brown   }
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* Create a sequential AIJ matrix on rank 0 convert it inplace to a new ViennaCL matrix, and apply it to a ViennaCL vector */
99c4762a1bSJed Brown   if (!rank) {
100c4762a1bSJed Brown     Mat               A;
101c4762a1bSJed Brown     Vec               v,r,v_vcl,r_vcl,d_vcl;
102c4762a1bSJed Brown     PetscInt          n=17,m=31,nz=5,i,cnt,j;
103c4762a1bSJed Brown     const PetscScalar val = 1.0;
104c4762a1bSJed Brown     PetscReal         dnorm;
105c4762a1bSJed Brown     const PetscReal   tol = 1e-5;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown     ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,nz,NULL,&A);CHKERRQ(ierr);
108c4762a1bSJed Brown 
109c4762a1bSJed Brown     /* Add nz arbitrary entries per row in arbitrary columns */
110c4762a1bSJed Brown     for (i=0;i<m;++i) {
111c4762a1bSJed Brown       for (cnt = 0; cnt<nz; ++cnt) {
112c4762a1bSJed Brown         j = (19 * cnt + (7*i + 3)) % n;
113c4762a1bSJed Brown         ierr = MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES);CHKERRQ(ierr);
114c4762a1bSJed Brown       }
115c4762a1bSJed Brown     }
116c4762a1bSJed Brown     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
117c4762a1bSJed Brown     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
118c4762a1bSJed Brown 
119c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix");CHKERRQ(ierr);
120c4762a1bSJed Brown 
121c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF,n,&v);CHKERRQ(ierr);
122c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF,m,&r);CHKERRQ(ierr);
123c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)r,"CPU result vector");CHKERRQ(ierr);
124c4762a1bSJed Brown     ierr = VecSet(v,val);CHKERRQ(ierr);
125c4762a1bSJed Brown     ierr = MatMult(A,v,r);CHKERRQ(ierr);
126c4762a1bSJed Brown 
127c4762a1bSJed Brown     ierr = MatConvert(A,MATSEQAIJVIENNACL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
128c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)A,"Converted ViennaCL Matrix");CHKERRQ(ierr);
129c4762a1bSJed Brown 
130c4762a1bSJed Brown     ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl);CHKERRQ(ierr);
131c4762a1bSJed Brown     ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl);CHKERRQ(ierr);
132c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector");CHKERRQ(ierr);
133c4762a1bSJed Brown     ierr = VecSet(v_vcl,val);CHKERRQ(ierr);
134c4762a1bSJed Brown     ierr = MatMult(A,v_vcl,r_vcl);CHKERRQ(ierr);
135c4762a1bSJed Brown 
136c4762a1bSJed Brown     ierr = VecDuplicate(r_vcl,&d_vcl);CHKERRQ(ierr);
137c4762a1bSJed Brown     ierr = VecCopy(r_vcl,d_vcl);CHKERRQ(ierr);
138*1e1ea65dSPierre Jolivet     ierr = VecAXPY(d_vcl,-1.0,r_vcl);CHKERRQ(ierr);
139c4762a1bSJed Brown     ierr = VecNorm(d_vcl,NORM_INFINITY,&dnorm);CHKERRQ(ierr);
140c4762a1bSJed Brown     if (dnorm > tol) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MPI CPU and MPI ViennaCL Vector results incompatible with inf norm greater than tolerance of %g",tol);
141c4762a1bSJed Brown 
142c4762a1bSJed Brown     ierr = VecDestroy(&v);CHKERRQ(ierr);
143c4762a1bSJed Brown     ierr = VecDestroy(&r);CHKERRQ(ierr);
144c4762a1bSJed Brown     ierr = VecDestroy(&v_vcl);CHKERRQ(ierr);
145c4762a1bSJed Brown     ierr = VecDestroy(&r_vcl);CHKERRQ(ierr);
146c4762a1bSJed Brown     ierr = VecDestroy(&d_vcl);CHKERRQ(ierr);
147c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
148c4762a1bSJed Brown   }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* Create a parallel AIJ matrix, convert it to a new ViennaCL matrix, and apply it to a ViennaCL vector */
151c4762a1bSJed Brown   if (size > 1) {
152c4762a1bSJed Brown     Mat               A,A_vcl;
153c4762a1bSJed Brown     Vec               v,r,v_vcl,r_vcl,d_vcl;
154c4762a1bSJed Brown     PetscInt          N=17,M=31,nz=5,i,cnt,j,rlow,rhigh;
155c4762a1bSJed Brown     const PetscScalar val = 1.0;
156c4762a1bSJed Brown     PetscReal         dnorm;
157c4762a1bSJed Brown     const PetscReal   tol=1e-5;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown     ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,M,N,nz,NULL,nz,NULL,&A);CHKERRQ(ierr);
160c4762a1bSJed Brown 
161c4762a1bSJed Brown     /* Add nz arbitrary entries per row in arbitrary columns */
162c4762a1bSJed Brown     ierr = MatGetOwnershipRange(A,&rlow,&rhigh);CHKERRQ(ierr);
163c4762a1bSJed Brown     for (i=rlow;i<rhigh;++i) {
164c4762a1bSJed Brown       for (cnt = 0; cnt<nz; ++cnt) {
165c4762a1bSJed Brown         j = (19 * cnt + (7*i + 3)) % N;
166*1e1ea65dSPierre Jolivet         ierr = MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES);CHKERRQ(ierr);
167c4762a1bSJed Brown       }
168c4762a1bSJed Brown     }
169c4762a1bSJed Brown     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
170c4762a1bSJed Brown     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
171c4762a1bSJed Brown 
172c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)A,"MPI CPU Matrix");CHKERRQ(ierr);
173c4762a1bSJed Brown 
174c4762a1bSJed Brown     ierr = MatCreateVecs(A,&v,&r);CHKERRQ(ierr);
175c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)r,"MPI CPU result vector");CHKERRQ(ierr);
176c4762a1bSJed Brown     ierr = VecSet(v,val);CHKERRQ(ierr);
177c4762a1bSJed Brown     ierr = MatMult(A,v,r);CHKERRQ(ierr);
178c4762a1bSJed Brown 
179c4762a1bSJed Brown     ierr = MatConvert(A,MATMPIAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl);CHKERRQ(ierr);
180c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)A_vcl,"New MPI ViennaCL Matrix");CHKERRQ(ierr);
181c4762a1bSJed Brown 
182c4762a1bSJed Brown     ierr = MatCreateVecs(A_vcl,&v_vcl,&r_vcl);CHKERRQ(ierr);
183c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector");CHKERRQ(ierr);
184c4762a1bSJed Brown     ierr = VecSet(v_vcl,val);CHKERRQ(ierr);
185c4762a1bSJed Brown     ierr = MatMult(A_vcl,v_vcl,r_vcl);CHKERRQ(ierr);
186c4762a1bSJed Brown 
187c4762a1bSJed Brown     ierr = VecDuplicate(r_vcl,&d_vcl);CHKERRQ(ierr);
188c4762a1bSJed Brown     ierr = VecCopy(r_vcl,d_vcl);CHKERRQ(ierr);
189*1e1ea65dSPierre Jolivet     ierr = VecAXPY(d_vcl,-1.0,r_vcl);CHKERRQ(ierr);
190c4762a1bSJed Brown     ierr = VecNorm(d_vcl,NORM_INFINITY,&dnorm);CHKERRQ(ierr);
191c4762a1bSJed Brown     if (dnorm > tol) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MPI CPU and MPI ViennaCL Vector results incompatible with inf norm greater than tolerance of %g",tol);
192c4762a1bSJed Brown 
193c4762a1bSJed Brown     ierr = VecDestroy(&v);CHKERRQ(ierr);
194c4762a1bSJed Brown     ierr = VecDestroy(&r);CHKERRQ(ierr);
195c4762a1bSJed Brown     ierr = VecDestroy(&v_vcl);CHKERRQ(ierr);
196c4762a1bSJed Brown     ierr = VecDestroy(&r_vcl);CHKERRQ(ierr);
197c4762a1bSJed Brown     ierr = VecDestroy(&d_vcl);CHKERRQ(ierr);
198c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
199c4762a1bSJed Brown     ierr = MatDestroy(&A_vcl);CHKERRQ(ierr);
200c4762a1bSJed Brown   }
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /* Create a parallel AIJ matrix, convert it in-place to a ViennaCL matrix, and apply it to a ViennaCL vector */
203c4762a1bSJed Brown   if (size > 1) {
204c4762a1bSJed Brown     Mat               A;
205c4762a1bSJed Brown     Vec               v,r,v_vcl,r_vcl,d_vcl;
206c4762a1bSJed Brown     PetscInt          N=17,M=31,nz=5,i,cnt,j,rlow,rhigh;
207c4762a1bSJed Brown     const PetscScalar val = 1.0;
208c4762a1bSJed Brown     PetscReal         dnorm;
209c4762a1bSJed Brown     const PetscReal   tol=1e-5;
210c4762a1bSJed Brown 
211c4762a1bSJed Brown     ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,M,N,nz,NULL,nz,NULL,&A);CHKERRQ(ierr);
212c4762a1bSJed Brown 
213c4762a1bSJed Brown     /* Add nz arbitrary entries per row in arbitrary columns */
214c4762a1bSJed Brown     ierr = MatGetOwnershipRange(A,&rlow,&rhigh);CHKERRQ(ierr);
215c4762a1bSJed Brown     for (i=rlow;i<rhigh;++i) {
216c4762a1bSJed Brown       for (cnt = 0; cnt<nz; ++cnt) {
217c4762a1bSJed Brown         j = (19 * cnt + (7*i + 3)) % N;
218c4762a1bSJed Brown         ierr = MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES);CHKERRQ(ierr);
219c4762a1bSJed Brown       }
220c4762a1bSJed Brown     }
221c4762a1bSJed Brown     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
222c4762a1bSJed Brown     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
223c4762a1bSJed Brown 
224c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)A,"MPI CPU Matrix");CHKERRQ(ierr);
225c4762a1bSJed Brown 
226c4762a1bSJed Brown     ierr = MatCreateVecs(A,&v,&r);CHKERRQ(ierr);
227c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)r,"MPI CPU result vector");CHKERRQ(ierr);
228c4762a1bSJed Brown     ierr = VecSet(v,val);CHKERRQ(ierr);
229c4762a1bSJed Brown     ierr = MatMult(A,v,r);CHKERRQ(ierr);
230c4762a1bSJed Brown 
231c4762a1bSJed Brown     ierr = MatConvert(A,MATMPIAIJVIENNACL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
232c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)A,"Converted MPI ViennaCL Matrix");CHKERRQ(ierr);
233c4762a1bSJed Brown 
234c4762a1bSJed Brown     ierr = MatCreateVecs(A,&v_vcl,&r_vcl);CHKERRQ(ierr);
235c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector");CHKERRQ(ierr);
236c4762a1bSJed Brown     ierr = VecSet(v_vcl,val);CHKERRQ(ierr);
237c4762a1bSJed Brown     ierr = MatMult(A,v_vcl,r_vcl);CHKERRQ(ierr);
238c4762a1bSJed Brown 
239c4762a1bSJed Brown     ierr = VecDuplicate(r_vcl,&d_vcl);CHKERRQ(ierr);
240c4762a1bSJed Brown     ierr = VecCopy(r_vcl,d_vcl);CHKERRQ(ierr);
241c4762a1bSJed Brown     ierr = VecAXPY(d_vcl,-1.0,r_vcl);CHKERRQ(ierr);
242c4762a1bSJed Brown     ierr = VecNorm(d_vcl,NORM_INFINITY,&dnorm);CHKERRQ(ierr);
243c4762a1bSJed Brown     if (dnorm > tol) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MPI CPU and MPI ViennaCL Vector results incompatible with inf norm greater than tolerance of %g",tol);
244c4762a1bSJed Brown 
245c4762a1bSJed Brown     ierr = VecDestroy(&v);CHKERRQ(ierr);
246c4762a1bSJed Brown     ierr = VecDestroy(&r);CHKERRQ(ierr);
247c4762a1bSJed Brown     ierr = VecDestroy(&v_vcl);CHKERRQ(ierr);
248c4762a1bSJed Brown     ierr = VecDestroy(&r_vcl);CHKERRQ(ierr);
249c4762a1bSJed Brown     ierr = VecDestroy(&d_vcl);CHKERRQ(ierr);
250c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
251c4762a1bSJed Brown   }
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   ierr = PetscFinalize();
254c4762a1bSJed Brown   return ierr;
255c4762a1bSJed Brown 
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown /*TEST
259c4762a1bSJed Brown 
260c4762a1bSJed Brown    build:
261c4762a1bSJed Brown       requires: viennacl
262c4762a1bSJed Brown 
263c4762a1bSJed Brown    test:
264c4762a1bSJed Brown       output_file: output/ex204.out
265c4762a1bSJed Brown 
266c4762a1bSJed Brown    test:
267c4762a1bSJed Brown       suffix: 2
268c4762a1bSJed Brown       nsize: 2
269c4762a1bSJed Brown       output_file: output/ex204.out
270c4762a1bSJed Brown 
271c4762a1bSJed Brown TEST*/
272