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 99566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 10c4762a1bSJed Brown 119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 13c4762a1bSJed Brown 14c4762a1bSJed Brown /* Construct a sequential ViennaCL matrix and vector */ 15dd400576SPatrick Sanan if (rank == 0) { 16c4762a1bSJed Brown Mat A_vcl; 17c4762a1bSJed Brown Vec v_vcl,r_vcl; 18c4762a1bSJed Brown PetscInt n = 17, m = 31,nz = 5,i,cnt,j; 19c4762a1bSJed Brown const PetscScalar val = 1.0; 20c4762a1bSJed Brown 219566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJViennaCL(PETSC_COMM_SELF,m,n,nz,NULL,&A_vcl)); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* Add nz arbitrary entries per row in arbitrary columns */ 24c4762a1bSJed Brown for (i=0;i<m;++i) { 25c4762a1bSJed Brown for (cnt = 0; cnt<nz; ++cnt) { 26c4762a1bSJed Brown j = (19 * cnt + (7*i + 3)) % n; 279566063dSJacob Faibussowitsch PetscCall(MatSetValue(A_vcl,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 28c4762a1bSJed Brown } 29c4762a1bSJed Brown } 309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A_vcl,MAT_FINAL_ASSEMBLY)); 319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A_vcl,MAT_FINAL_ASSEMBLY)); 32c4762a1bSJed Brown 339566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 349566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 359566063dSJacob Faibussowitsch PetscCall(VecSet(v_vcl,val)); 36c4762a1bSJed Brown 379566063dSJacob Faibussowitsch PetscCall(MatMult(A_vcl,v_vcl,r_vcl)); 38c4762a1bSJed Brown 399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_vcl)); 409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r_vcl)); 419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_vcl)); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* Create a sequential AIJ matrix on rank 0 convert it to a new ViennaCL matrix, and apply it to a ViennaCL vector */ 45dd400576SPatrick Sanan if (rank == 0) { 46c4762a1bSJed Brown Mat A,A_vcl; 47c4762a1bSJed Brown Vec v,r,v_vcl,r_vcl,d_vcl; 48c4762a1bSJed Brown PetscInt n=17,m=31,nz=5,i,cnt,j; 49c4762a1bSJed Brown const PetscScalar val = 1.0; 50c4762a1bSJed Brown PetscReal dnorm; 51c4762a1bSJed Brown const PetscReal tol = 1e-5; 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,nz,NULL,&A)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Add nz arbitrary entries per row in arbitrary columns */ 56c4762a1bSJed Brown for (i=0;i<m;++i) { 57c4762a1bSJed Brown for (cnt = 0; cnt<nz; ++cnt) { 58c4762a1bSJed Brown j = (19 * cnt + (7*i + 3)) % n; 599566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,j,(PetscScalar) (0.3 * i + j),INSERT_VALUES)); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown } 629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix")); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&v)); 689566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&r)); 699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r,"CPU result vector")); 709566063dSJacob Faibussowitsch PetscCall(VecSet(v,val)); 719566063dSJacob Faibussowitsch PetscCall(MatMult(A,v,r)); 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATSEQAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl)); 749566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A_vcl,"New ViennaCL Matrix")); 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 779566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 789566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 799566063dSJacob Faibussowitsch PetscCall(VecSet(v_vcl,val)); 809566063dSJacob Faibussowitsch PetscCall(MatMult(A_vcl,v_vcl,r_vcl)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r_vcl,&d_vcl)); 839566063dSJacob Faibussowitsch PetscCall(VecCopy(r_vcl,d_vcl)); 849566063dSJacob Faibussowitsch PetscCall(VecAXPY(d_vcl,-1.0,r_vcl)); 859566063dSJacob Faibussowitsch PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 86*08401ef6SPierre 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); 87c4762a1bSJed Brown 889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_vcl)); 919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r_vcl)); 929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d_vcl)); 939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_vcl)); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Create a sequential AIJ matrix on rank 0 convert it inplace to a new ViennaCL matrix, and apply it to a ViennaCL vector */ 98dd400576SPatrick Sanan if (rank == 0) { 99c4762a1bSJed Brown Mat A; 100c4762a1bSJed Brown Vec v,r,v_vcl,r_vcl,d_vcl; 101c4762a1bSJed Brown PetscInt n=17,m=31,nz=5,i,cnt,j; 102c4762a1bSJed Brown const PetscScalar val = 1.0; 103c4762a1bSJed Brown PetscReal dnorm; 104c4762a1bSJed Brown const PetscReal tol = 1e-5; 105c4762a1bSJed Brown 1069566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,nz,NULL,&A)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Add nz arbitrary entries per row in arbitrary columns */ 109c4762a1bSJed Brown for (i=0;i<m;++i) { 110c4762a1bSJed Brown for (cnt = 0; cnt<nz; ++cnt) { 111c4762a1bSJed Brown j = (19 * cnt + (7*i + 3)) % n; 1129566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown } 1159566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1169566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 117c4762a1bSJed Brown 1189566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix")); 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&v)); 1219566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&r)); 1229566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r,"CPU result vector")); 1239566063dSJacob Faibussowitsch PetscCall(VecSet(v,val)); 1249566063dSJacob Faibussowitsch PetscCall(MatMult(A,v,r)); 125c4762a1bSJed Brown 1269566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATSEQAIJVIENNACL,MAT_INPLACE_MATRIX,&A)); 1279566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"Converted ViennaCL Matrix")); 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 1309566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 1319566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 1329566063dSJacob Faibussowitsch PetscCall(VecSet(v_vcl,val)); 1339566063dSJacob Faibussowitsch PetscCall(MatMult(A,v_vcl,r_vcl)); 134c4762a1bSJed Brown 1359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r_vcl,&d_vcl)); 1369566063dSJacob Faibussowitsch PetscCall(VecCopy(r_vcl,d_vcl)); 1379566063dSJacob Faibussowitsch PetscCall(VecAXPY(d_vcl,-1.0,r_vcl)); 1389566063dSJacob Faibussowitsch PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 139*08401ef6SPierre 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); 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_vcl)); 1449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r_vcl)); 1459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d_vcl)); 1469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Create a parallel AIJ matrix, convert it to a new ViennaCL matrix, and apply it to a ViennaCL vector */ 150c4762a1bSJed Brown if (size > 1) { 151c4762a1bSJed Brown Mat A,A_vcl; 152c4762a1bSJed Brown Vec v,r,v_vcl,r_vcl,d_vcl; 153c4762a1bSJed Brown PetscInt N=17,M=31,nz=5,i,cnt,j,rlow,rhigh; 154c4762a1bSJed Brown const PetscScalar val = 1.0; 155c4762a1bSJed Brown PetscReal dnorm; 156c4762a1bSJed Brown const PetscReal tol=1e-5; 157c4762a1bSJed Brown 1589566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,M,N,nz,NULL,nz,NULL,&A)); 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* Add nz arbitrary entries per row in arbitrary columns */ 1619566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rlow,&rhigh)); 162c4762a1bSJed Brown for (i=rlow;i<rhigh;++i) { 163c4762a1bSJed Brown for (cnt = 0; cnt<nz; ++cnt) { 164c4762a1bSJed Brown j = (19 * cnt + (7*i + 3)) % N; 1659566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown } 1689566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1699566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 170c4762a1bSJed Brown 1719566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"MPI CPU Matrix")); 172c4762a1bSJed Brown 1739566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&v,&r)); 1749566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r,"MPI CPU result vector")); 1759566063dSJacob Faibussowitsch PetscCall(VecSet(v,val)); 1769566063dSJacob Faibussowitsch PetscCall(MatMult(A,v,r)); 177c4762a1bSJed Brown 1789566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATMPIAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl)); 1799566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A_vcl,"New MPI ViennaCL Matrix")); 180c4762a1bSJed Brown 1819566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A_vcl,&v_vcl,&r_vcl)); 1829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 1839566063dSJacob Faibussowitsch PetscCall(VecSet(v_vcl,val)); 1849566063dSJacob Faibussowitsch PetscCall(MatMult(A_vcl,v_vcl,r_vcl)); 185c4762a1bSJed Brown 1869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r_vcl,&d_vcl)); 1879566063dSJacob Faibussowitsch PetscCall(VecCopy(r_vcl,d_vcl)); 1889566063dSJacob Faibussowitsch PetscCall(VecAXPY(d_vcl,-1.0,r_vcl)); 1899566063dSJacob Faibussowitsch PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 190*08401ef6SPierre 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); 191c4762a1bSJed Brown 1929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_vcl)); 1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r_vcl)); 1969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d_vcl)); 1979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_vcl)); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* Create a parallel AIJ matrix, convert it in-place to a ViennaCL matrix, and apply it to a ViennaCL vector */ 202c4762a1bSJed Brown if (size > 1) { 203c4762a1bSJed Brown Mat A; 204c4762a1bSJed Brown Vec v,r,v_vcl,r_vcl,d_vcl; 205c4762a1bSJed Brown PetscInt N=17,M=31,nz=5,i,cnt,j,rlow,rhigh; 206c4762a1bSJed Brown const PetscScalar val = 1.0; 207c4762a1bSJed Brown PetscReal dnorm; 208c4762a1bSJed Brown const PetscReal tol=1e-5; 209c4762a1bSJed Brown 2109566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,M,N,nz,NULL,nz,NULL,&A)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* Add nz arbitrary entries per row in arbitrary columns */ 2139566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rlow,&rhigh)); 214c4762a1bSJed Brown for (i=rlow;i<rhigh;++i) { 215c4762a1bSJed Brown for (cnt = 0; cnt<nz; ++cnt) { 216c4762a1bSJed Brown j = (19 * cnt + (7*i + 3)) % N; 2179566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown } 2209566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 2219566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 222c4762a1bSJed Brown 2239566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"MPI CPU Matrix")); 224c4762a1bSJed Brown 2259566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&v,&r)); 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r,"MPI CPU result vector")); 2279566063dSJacob Faibussowitsch PetscCall(VecSet(v,val)); 2289566063dSJacob Faibussowitsch PetscCall(MatMult(A,v,r)); 229c4762a1bSJed Brown 2309566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATMPIAIJVIENNACL,MAT_INPLACE_MATRIX,&A)); 2319566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"Converted MPI ViennaCL Matrix")); 232c4762a1bSJed Brown 2339566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&v_vcl,&r_vcl)); 2349566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 2359566063dSJacob Faibussowitsch PetscCall(VecSet(v_vcl,val)); 2369566063dSJacob Faibussowitsch PetscCall(MatMult(A,v_vcl,r_vcl)); 237c4762a1bSJed Brown 2389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r_vcl,&d_vcl)); 2399566063dSJacob Faibussowitsch PetscCall(VecCopy(r_vcl,d_vcl)); 2409566063dSJacob Faibussowitsch PetscCall(VecAXPY(d_vcl,-1.0,r_vcl)); 2419566063dSJacob Faibussowitsch PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 242*08401ef6SPierre 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); 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 2459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 2469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_vcl)); 2479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r_vcl)); 2489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d_vcl)); 2499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 2529566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 253b122ec5aSJacob Faibussowitsch return 0; 254c4762a1bSJed Brown 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown /*TEST 258c4762a1bSJed Brown 259c4762a1bSJed Brown build: 260c4762a1bSJed Brown requires: viennacl 261c4762a1bSJed Brown 262c4762a1bSJed Brown test: 263c4762a1bSJed Brown output_file: output/ex204.out 264c4762a1bSJed Brown 265c4762a1bSJed Brown test: 266c4762a1bSJed Brown suffix: 2 267c4762a1bSJed Brown nsize: 2 268c4762a1bSJed Brown output_file: output/ex204.out 269c4762a1bSJed Brown 270c4762a1bSJed Brown TEST*/ 271