xref: /petsc/src/mat/tests/ex204.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscMPIInt rank,size;
8c4762a1bSJed Brown 
9*327415f7SBarry Smith   PetscFunctionBeginUser;
109566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
11c4762a1bSJed Brown 
129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   /* Construct a sequential ViennaCL matrix and vector */
16dd400576SPatrick Sanan   if (rank == 0) {
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 
229566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJViennaCL(PETSC_COMM_SELF,m,n,nz,NULL,&A_vcl));
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;
289566063dSJacob Faibussowitsch         PetscCall(MatSetValue(A_vcl,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES));
29c4762a1bSJed Brown       }
30c4762a1bSJed Brown     }
319566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A_vcl,MAT_FINAL_ASSEMBLY));
329566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A_vcl,MAT_FINAL_ASSEMBLY));
33c4762a1bSJed Brown 
349566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl));
359566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl));
369566063dSJacob Faibussowitsch     PetscCall(VecSet(v_vcl,val));
37c4762a1bSJed Brown 
389566063dSJacob Faibussowitsch     PetscCall(MatMult(A_vcl,v_vcl,r_vcl));
39c4762a1bSJed Brown 
409566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v_vcl));
419566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r_vcl));
429566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A_vcl));
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 */
46dd400576SPatrick Sanan   if (rank == 0) {
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 
549566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,nz,NULL,&A));
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;
609566063dSJacob Faibussowitsch         PetscCall(MatSetValue(A,i,j,(PetscScalar) (0.3 * i + j),INSERT_VALUES));
61c4762a1bSJed Brown       }
62c4762a1bSJed Brown     }
639566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
649566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix"));
67c4762a1bSJed Brown 
689566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&v));
699566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&r));
709566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r,"CPU result vector"));
719566063dSJacob Faibussowitsch     PetscCall(VecSet(v,val));
729566063dSJacob Faibussowitsch     PetscCall(MatMult(A,v,r));
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATSEQAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl));
759566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)A_vcl,"New ViennaCL Matrix"));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl));
789566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl));
799566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector"));
809566063dSJacob Faibussowitsch     PetscCall(VecSet(v_vcl,val));
819566063dSJacob Faibussowitsch     PetscCall(MatMult(A_vcl,v_vcl,r_vcl));
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(r_vcl,&d_vcl));
849566063dSJacob Faibussowitsch     PetscCall(VecCopy(r_vcl,d_vcl));
859566063dSJacob Faibussowitsch     PetscCall(VecAXPY(d_vcl,-1.0,r_vcl));
869566063dSJacob Faibussowitsch     PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm));
8708401ef6SPierre Jolivet     PetscCheck(dnorm <= tol,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 
899566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v));
909566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
919566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v_vcl));
929566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r_vcl));
939566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&d_vcl));
949566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
959566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A_vcl));
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 */
99dd400576SPatrick Sanan   if (rank == 0) {
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 
1079566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,nz,NULL,&A));
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;
1139566063dSJacob Faibussowitsch         PetscCall(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES));
114c4762a1bSJed Brown       }
115c4762a1bSJed Brown     }
1169566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1179566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
118c4762a1bSJed Brown 
1199566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix"));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&v));
1229566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&r));
1239566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r,"CPU result vector"));
1249566063dSJacob Faibussowitsch     PetscCall(VecSet(v,val));
1259566063dSJacob Faibussowitsch     PetscCall(MatMult(A,v,r));
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATSEQAIJVIENNACL,MAT_INPLACE_MATRIX,&A));
1289566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)A,"Converted ViennaCL Matrix"));
129c4762a1bSJed Brown 
1309566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl));
1319566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl));
1329566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector"));
1339566063dSJacob Faibussowitsch     PetscCall(VecSet(v_vcl,val));
1349566063dSJacob Faibussowitsch     PetscCall(MatMult(A,v_vcl,r_vcl));
135c4762a1bSJed Brown 
1369566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(r_vcl,&d_vcl));
1379566063dSJacob Faibussowitsch     PetscCall(VecCopy(r_vcl,d_vcl));
1389566063dSJacob Faibussowitsch     PetscCall(VecAXPY(d_vcl,-1.0,r_vcl));
1399566063dSJacob Faibussowitsch     PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm));
14008401ef6SPierre Jolivet     PetscCheck(dnorm <= tol,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 
1429566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v));
1439566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
1449566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v_vcl));
1459566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r_vcl));
1469566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&d_vcl));
1479566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
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 
1599566063dSJacob Faibussowitsch     PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,M,N,nz,NULL,nz,NULL,&A));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown     /* Add nz arbitrary entries per row in arbitrary columns */
1629566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A,&rlow,&rhigh));
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;
1669566063dSJacob Faibussowitsch         PetscCall(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES));
167c4762a1bSJed Brown       }
168c4762a1bSJed Brown     }
1699566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1709566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)A,"MPI CPU Matrix"));
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A,&v,&r));
1759566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r,"MPI CPU result vector"));
1769566063dSJacob Faibussowitsch     PetscCall(VecSet(v,val));
1779566063dSJacob Faibussowitsch     PetscCall(MatMult(A,v,r));
178c4762a1bSJed Brown 
1799566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATMPIAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl));
1809566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)A_vcl,"New MPI ViennaCL Matrix"));
181c4762a1bSJed Brown 
1829566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A_vcl,&v_vcl,&r_vcl));
1839566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector"));
1849566063dSJacob Faibussowitsch     PetscCall(VecSet(v_vcl,val));
1859566063dSJacob Faibussowitsch     PetscCall(MatMult(A_vcl,v_vcl,r_vcl));
186c4762a1bSJed Brown 
1879566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(r_vcl,&d_vcl));
1889566063dSJacob Faibussowitsch     PetscCall(VecCopy(r_vcl,d_vcl));
1899566063dSJacob Faibussowitsch     PetscCall(VecAXPY(d_vcl,-1.0,r_vcl));
1909566063dSJacob Faibussowitsch     PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm));
19108401ef6SPierre Jolivet     PetscCheck(dnorm <= tol,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 
1939566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v));
1949566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
1959566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v_vcl));
1969566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r_vcl));
1979566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&d_vcl));
1989566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
1999566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A_vcl));
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 
2119566063dSJacob Faibussowitsch     PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,M,N,nz,NULL,nz,NULL,&A));
212c4762a1bSJed Brown 
213c4762a1bSJed Brown     /* Add nz arbitrary entries per row in arbitrary columns */
2149566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A,&rlow,&rhigh));
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;
2189566063dSJacob Faibussowitsch         PetscCall(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES));
219c4762a1bSJed Brown       }
220c4762a1bSJed Brown     }
2219566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
2229566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
223c4762a1bSJed Brown 
2249566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)A,"MPI CPU Matrix"));
225c4762a1bSJed Brown 
2269566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A,&v,&r));
2279566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r,"MPI CPU result vector"));
2289566063dSJacob Faibussowitsch     PetscCall(VecSet(v,val));
2299566063dSJacob Faibussowitsch     PetscCall(MatMult(A,v,r));
230c4762a1bSJed Brown 
2319566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATMPIAIJVIENNACL,MAT_INPLACE_MATRIX,&A));
2329566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)A,"Converted MPI ViennaCL Matrix"));
233c4762a1bSJed Brown 
2349566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A,&v_vcl,&r_vcl));
2359566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector"));
2369566063dSJacob Faibussowitsch     PetscCall(VecSet(v_vcl,val));
2379566063dSJacob Faibussowitsch     PetscCall(MatMult(A,v_vcl,r_vcl));
238c4762a1bSJed Brown 
2399566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(r_vcl,&d_vcl));
2409566063dSJacob Faibussowitsch     PetscCall(VecCopy(r_vcl,d_vcl));
2419566063dSJacob Faibussowitsch     PetscCall(VecAXPY(d_vcl,-1.0,r_vcl));
2429566063dSJacob Faibussowitsch     PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm));
24308401ef6SPierre Jolivet     PetscCheck(dnorm <= tol,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 
2459566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v));
2469566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
2479566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v_vcl));
2489566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r_vcl));
2499566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&d_vcl));
2509566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
251c4762a1bSJed Brown   }
252c4762a1bSJed Brown 
2539566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
254b122ec5aSJacob Faibussowitsch   return 0;
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