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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 10c4762a1bSJed Brown 11*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 12*9566063dSJacob 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 21*9566063dSJacob 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; 27*9566063dSJacob Faibussowitsch PetscCall(MatSetValue(A_vcl,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 28c4762a1bSJed Brown } 29c4762a1bSJed Brown } 30*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A_vcl,MAT_FINAL_ASSEMBLY)); 31*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A_vcl,MAT_FINAL_ASSEMBLY)); 32c4762a1bSJed Brown 33*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 34*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 35*9566063dSJacob Faibussowitsch PetscCall(VecSet(v_vcl,val)); 36c4762a1bSJed Brown 37*9566063dSJacob Faibussowitsch PetscCall(MatMult(A_vcl,v_vcl,r_vcl)); 38c4762a1bSJed Brown 39*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_vcl)); 40*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r_vcl)); 41*9566063dSJacob 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 53*9566063dSJacob 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; 59*9566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,j,(PetscScalar) (0.3 * i + j),INSERT_VALUES)); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown } 62*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 63*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 64c4762a1bSJed Brown 65*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix")); 66c4762a1bSJed Brown 67*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&v)); 68*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&r)); 69*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r,"CPU result vector")); 70*9566063dSJacob Faibussowitsch PetscCall(VecSet(v,val)); 71*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,v,r)); 72c4762a1bSJed Brown 73*9566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATSEQAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl)); 74*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A_vcl,"New ViennaCL Matrix")); 75c4762a1bSJed Brown 76*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 77*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 78*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 79*9566063dSJacob Faibussowitsch PetscCall(VecSet(v_vcl,val)); 80*9566063dSJacob Faibussowitsch PetscCall(MatMult(A_vcl,v_vcl,r_vcl)); 81c4762a1bSJed Brown 82*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r_vcl,&d_vcl)); 83*9566063dSJacob Faibussowitsch PetscCall(VecCopy(r_vcl,d_vcl)); 84*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(d_vcl,-1.0,r_vcl)); 85*9566063dSJacob Faibussowitsch PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 862c71b3e2SJacob Faibussowitsch PetscCheckFalse(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 88*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 89*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 90*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_vcl)); 91*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r_vcl)); 92*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d_vcl)); 93*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 94*9566063dSJacob 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 106*9566063dSJacob 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; 112*9566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown } 115*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 116*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 117c4762a1bSJed Brown 118*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix")); 119c4762a1bSJed Brown 120*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&v)); 121*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&r)); 122*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r,"CPU result vector")); 123*9566063dSJacob Faibussowitsch PetscCall(VecSet(v,val)); 124*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,v,r)); 125c4762a1bSJed Brown 126*9566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATSEQAIJVIENNACL,MAT_INPLACE_MATRIX,&A)); 127*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"Converted ViennaCL Matrix")); 128c4762a1bSJed Brown 129*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 130*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 131*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 132*9566063dSJacob Faibussowitsch PetscCall(VecSet(v_vcl,val)); 133*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,v_vcl,r_vcl)); 134c4762a1bSJed Brown 135*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r_vcl,&d_vcl)); 136*9566063dSJacob Faibussowitsch PetscCall(VecCopy(r_vcl,d_vcl)); 137*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(d_vcl,-1.0,r_vcl)); 138*9566063dSJacob Faibussowitsch PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 1392c71b3e2SJacob Faibussowitsch PetscCheckFalse(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 141*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 142*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 143*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_vcl)); 144*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r_vcl)); 145*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d_vcl)); 146*9566063dSJacob 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 158*9566063dSJacob 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 */ 161*9566063dSJacob 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; 165*9566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown } 168*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 169*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 170c4762a1bSJed Brown 171*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"MPI CPU Matrix")); 172c4762a1bSJed Brown 173*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&v,&r)); 174*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r,"MPI CPU result vector")); 175*9566063dSJacob Faibussowitsch PetscCall(VecSet(v,val)); 176*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,v,r)); 177c4762a1bSJed Brown 178*9566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATMPIAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl)); 179*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A_vcl,"New MPI ViennaCL Matrix")); 180c4762a1bSJed Brown 181*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A_vcl,&v_vcl,&r_vcl)); 182*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 183*9566063dSJacob Faibussowitsch PetscCall(VecSet(v_vcl,val)); 184*9566063dSJacob Faibussowitsch PetscCall(MatMult(A_vcl,v_vcl,r_vcl)); 185c4762a1bSJed Brown 186*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r_vcl,&d_vcl)); 187*9566063dSJacob Faibussowitsch PetscCall(VecCopy(r_vcl,d_vcl)); 188*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(d_vcl,-1.0,r_vcl)); 189*9566063dSJacob Faibussowitsch PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 1902c71b3e2SJacob Faibussowitsch PetscCheckFalse(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 192*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 193*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 194*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_vcl)); 195*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r_vcl)); 196*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d_vcl)); 197*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 198*9566063dSJacob 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 210*9566063dSJacob 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 */ 213*9566063dSJacob 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; 217*9566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown } 220*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 221*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 222c4762a1bSJed Brown 223*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"MPI CPU Matrix")); 224c4762a1bSJed Brown 225*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&v,&r)); 226*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r,"MPI CPU result vector")); 227*9566063dSJacob Faibussowitsch PetscCall(VecSet(v,val)); 228*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,v,r)); 229c4762a1bSJed Brown 230*9566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATMPIAIJVIENNACL,MAT_INPLACE_MATRIX,&A)); 231*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"Converted MPI ViennaCL Matrix")); 232c4762a1bSJed Brown 233*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&v_vcl,&r_vcl)); 234*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 235*9566063dSJacob Faibussowitsch PetscCall(VecSet(v_vcl,val)); 236*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,v_vcl,r_vcl)); 237c4762a1bSJed Brown 238*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r_vcl,&d_vcl)); 239*9566063dSJacob Faibussowitsch PetscCall(VecCopy(r_vcl,d_vcl)); 240*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(d_vcl,-1.0,r_vcl)); 241*9566063dSJacob Faibussowitsch PetscCall(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 2422c71b3e2SJacob Faibussowitsch PetscCheckFalse(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 244*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 245*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 246*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_vcl)); 247*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r_vcl)); 248*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d_vcl)); 249*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252*9566063dSJacob 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