static char help[] = "Test ViennaCL Matrix Conversions"; #include int main(int argc,char **args) { PetscMPIInt rank,size; PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); /* Construct a sequential ViennaCL matrix and vector */ if (rank == 0) { Mat A_vcl; Vec v_vcl,r_vcl; PetscInt n = 17, m = 31,nz = 5,i,cnt,j; const PetscScalar val = 1.0; PetscCall(MatCreateSeqAIJViennaCL(PETSC_COMM_SELF,m,n,nz,NULL,&A_vcl)); /* Add nz arbitrary entries per row in arbitrary columns */ for (i=0;i 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); PetscCall(VecDestroy(&v)); PetscCall(VecDestroy(&r)); PetscCall(VecDestroy(&v_vcl)); PetscCall(VecDestroy(&r_vcl)); PetscCall(VecDestroy(&d_vcl)); PetscCall(MatDestroy(&A)); PetscCall(MatDestroy(&A_vcl)); } /* Create a sequential AIJ matrix on rank 0 convert it inplace to a new ViennaCL matrix, and apply it to a ViennaCL vector */ if (rank == 0) { Mat A; Vec v,r,v_vcl,r_vcl,d_vcl; PetscInt n=17,m=31,nz=5,i,cnt,j; const PetscScalar val = 1.0; PetscReal dnorm; const PetscReal tol = 1e-5; PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,nz,NULL,&A)); /* Add nz arbitrary entries per row in arbitrary columns */ for (i=0;i 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); PetscCall(VecDestroy(&v)); PetscCall(VecDestroy(&r)); PetscCall(VecDestroy(&v_vcl)); PetscCall(VecDestroy(&r_vcl)); PetscCall(VecDestroy(&d_vcl)); PetscCall(MatDestroy(&A)); } /* Create a parallel AIJ matrix, convert it to a new ViennaCL matrix, and apply it to a ViennaCL vector */ if (size > 1) { Mat A,A_vcl; Vec v,r,v_vcl,r_vcl,d_vcl; PetscInt N=17,M=31,nz=5,i,cnt,j,rlow,rhigh; const PetscScalar val = 1.0; PetscReal dnorm; const PetscReal tol=1e-5; PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,M,N,nz,NULL,nz,NULL,&A)); /* Add nz arbitrary entries per row in arbitrary columns */ PetscCall(MatGetOwnershipRange(A,&rlow,&rhigh)); for (i=rlow;i 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); PetscCall(VecDestroy(&v)); PetscCall(VecDestroy(&r)); PetscCall(VecDestroy(&v_vcl)); PetscCall(VecDestroy(&r_vcl)); PetscCall(VecDestroy(&d_vcl)); PetscCall(MatDestroy(&A)); PetscCall(MatDestroy(&A_vcl)); } /* Create a parallel AIJ matrix, convert it in-place to a ViennaCL matrix, and apply it to a ViennaCL vector */ if (size > 1) { Mat A; Vec v,r,v_vcl,r_vcl,d_vcl; PetscInt N=17,M=31,nz=5,i,cnt,j,rlow,rhigh; const PetscScalar val = 1.0; PetscReal dnorm; const PetscReal tol=1e-5; PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,M,N,nz,NULL,nz,NULL,&A)); /* Add nz arbitrary entries per row in arbitrary columns */ PetscCall(MatGetOwnershipRange(A,&rlow,&rhigh)); for (i=rlow;i 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); PetscCall(VecDestroy(&v)); PetscCall(VecDestroy(&r)); PetscCall(VecDestroy(&v_vcl)); PetscCall(VecDestroy(&r_vcl)); PetscCall(VecDestroy(&d_vcl)); PetscCall(MatDestroy(&A)); } PetscCall(PetscFinalize()); return 0; } /*TEST build: requires: viennacl test: output_file: output/ex204.out test: suffix: 2 nsize: 2 output_file: output/ex204.out TEST*/