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 12*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 13*5f80ce2aSJacob Faibussowitsch CHKERRMPI(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 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A_vcl,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown } 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A_vcl,MAT_FINAL_ASSEMBLY)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A_vcl,MAT_FINAL_ASSEMBLY)); 33c4762a1bSJed Brown 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_vcl,val)); 37c4762a1bSJed Brown 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A_vcl,v_vcl,r_vcl)); 39c4762a1bSJed Brown 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_vcl)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r_vcl)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,i,j,(PetscScalar) (0.3 * i + j),INSERT_VALUES)); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown } 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 65c4762a1bSJed Brown 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix")); 67c4762a1bSJed Brown 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&v)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&r)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r,"CPU result vector")); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,val)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v,r)); 73c4762a1bSJed Brown 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSEQAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A_vcl,"New ViennaCL Matrix")); 76c4762a1bSJed Brown 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_vcl,val)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A_vcl,v_vcl,r_vcl)); 82c4762a1bSJed Brown 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(r_vcl,&d_vcl)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(r_vcl,d_vcl)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(d_vcl,-1.0,r_vcl)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 872c71b3e2SJacob 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); 88c4762a1bSJed Brown 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_vcl)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r_vcl)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d_vcl)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown } 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 118c4762a1bSJed Brown 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix")); 120c4762a1bSJed Brown 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&v)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&r)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r,"CPU result vector")); 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,val)); 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v,r)); 126c4762a1bSJed Brown 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSEQAIJVIENNACL,MAT_INPLACE_MATRIX,&A)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"Converted ViennaCL Matrix")); 129c4762a1bSJed Brown 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_vcl,val)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v_vcl,r_vcl)); 135c4762a1bSJed Brown 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(r_vcl,&d_vcl)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(r_vcl,d_vcl)); 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(d_vcl,-1.0,r_vcl)); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 1402c71b3e2SJacob 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); 141c4762a1bSJed Brown 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_vcl)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r_vcl)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d_vcl)); 147*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown } 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 171c4762a1bSJed Brown 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"MPI CPU Matrix")); 173c4762a1bSJed Brown 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&v,&r)); 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r,"MPI CPU result vector")); 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,val)); 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v,r)); 178c4762a1bSJed Brown 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATMPIAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl)); 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A_vcl,"New MPI ViennaCL Matrix")); 181c4762a1bSJed Brown 182*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A_vcl,&v_vcl,&r_vcl)); 183*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 184*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_vcl,val)); 185*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A_vcl,v_vcl,r_vcl)); 186c4762a1bSJed Brown 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(r_vcl,&d_vcl)); 188*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(r_vcl,d_vcl)); 189*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(d_vcl,-1.0,r_vcl)); 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 1912c71b3e2SJacob 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); 192c4762a1bSJed Brown 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 194*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_vcl)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r_vcl)); 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d_vcl)); 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 199*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 214*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown } 221*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 222*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 223c4762a1bSJed Brown 224*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"MPI CPU Matrix")); 225c4762a1bSJed Brown 226*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&v,&r)); 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r,"MPI CPU result vector")); 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,val)); 229*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v,r)); 230c4762a1bSJed Brown 231*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATMPIAIJVIENNACL,MAT_INPLACE_MATRIX,&A)); 232*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"Converted MPI ViennaCL Matrix")); 233c4762a1bSJed Brown 234*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&v_vcl,&r_vcl)); 235*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 236*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_vcl,val)); 237*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v_vcl,r_vcl)); 238c4762a1bSJed Brown 239*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(r_vcl,&d_vcl)); 240*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(r_vcl,d_vcl)); 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(d_vcl,-1.0,r_vcl)); 242*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(d_vcl,NORM_INFINITY,&dnorm)); 2432c71b3e2SJacob 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); 244c4762a1bSJed Brown 245*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 247*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_vcl)); 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r_vcl)); 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d_vcl)); 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 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