1*c4762a1bSJed Brown static char help[] = "Test ViennaCL Matrix Conversions"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown int main(int argc,char **args) 6*c4762a1bSJed Brown { 7*c4762a1bSJed Brown PetscErrorCode ierr; 8*c4762a1bSJed Brown PetscMPIInt rank,size; 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 13*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown /* Construct a sequential ViennaCL matrix and vector */ 16*c4762a1bSJed Brown if (!rank) { 17*c4762a1bSJed Brown Mat A_vcl; 18*c4762a1bSJed Brown Vec v_vcl,r_vcl; 19*c4762a1bSJed Brown PetscInt n = 17, m = 31,nz = 5,i,cnt,j; 20*c4762a1bSJed Brown const PetscScalar val = 1.0; 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown ierr = MatCreateSeqAIJViennaCL(PETSC_COMM_SELF,m,n,nz,NULL,&A_vcl);CHKERRQ(ierr); 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown /* Add nz arbitrary entries per row in arbitrary columns */ 25*c4762a1bSJed Brown for(i=0;i<m;++i){ 26*c4762a1bSJed Brown for(cnt = 0; cnt<nz; ++cnt) { 27*c4762a1bSJed Brown j = (19 * cnt + (7*i + 3)) % n; 28*c4762a1bSJed Brown ierr = MatSetValue(A_vcl,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES);CHKERRQ(ierr); 29*c4762a1bSJed Brown } 30*c4762a1bSJed Brown } 31*c4762a1bSJed Brown ierr = MatAssemblyBegin(A_vcl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = MatAssemblyEnd(A_vcl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = VecSet(v_vcl,val);CHKERRQ(ierr); 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown ierr = MatMult(A_vcl,v_vcl,r_vcl);CHKERRQ(ierr); 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown ierr = VecDestroy(&v_vcl);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = VecDestroy(&r_vcl);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = MatDestroy(&A_vcl);CHKERRQ(ierr); 43*c4762a1bSJed Brown } 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown /* Create a sequential AIJ matrix on rank 0 convert it to a new ViennaCL matrix, and apply it to a ViennaCL vector */ 46*c4762a1bSJed Brown if(!rank) { 47*c4762a1bSJed Brown Mat A,A_vcl; 48*c4762a1bSJed Brown Vec v,r,v_vcl,r_vcl,d_vcl; 49*c4762a1bSJed Brown PetscInt n=17,m=31,nz=5,i,cnt,j; 50*c4762a1bSJed Brown const PetscScalar val = 1.0; 51*c4762a1bSJed Brown PetscReal dnorm; 52*c4762a1bSJed Brown const PetscReal tol = 1e-5; 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,nz,NULL,&A);CHKERRQ(ierr); 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown /* Add nz arbitrary entries per row in arbitrary columns */ 57*c4762a1bSJed Brown for(i=0;i<m;++i){ 58*c4762a1bSJed Brown for(cnt = 0; cnt<nz; ++cnt){ 59*c4762a1bSJed Brown j = (19 * cnt + (7*i + 3)) % n; 60*c4762a1bSJed Brown ierr = MatSetValue(A,i,j,(PetscScalar) (0.3 * i + j),INSERT_VALUES);CHKERRQ(ierr); 61*c4762a1bSJed Brown } 62*c4762a1bSJed Brown } 63*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix");CHKERRQ(ierr); 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,n,&v);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,m,&r);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)r,"CPU result vector");CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = VecSet(v,val);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = MatMult(A,v,r);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = MatConvert(A,MATSEQAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)A_vcl,"New ViennaCL Matrix");CHKERRQ(ierr); 76*c4762a1bSJed Brown 77*c4762a1bSJed Brown ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector");CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = VecSet(v_vcl,val);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = MatMult(A_vcl,v_vcl,r_vcl);CHKERRQ(ierr); 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown ierr = VecDuplicate(r_vcl,&d_vcl);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = VecCopy(r_vcl,d_vcl);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = VecAXPY(d_vcl,-1.0,r_vcl); 86*c4762a1bSJed Brown ierr = VecNorm(d_vcl,NORM_INFINITY,&dnorm);CHKERRQ(ierr); 87*c4762a1bSJed Brown if (dnorm > tol) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"Sequential CPU and MPI ViennaCL vector results incompatible with inf norm greater than tolerance of %g",tol); 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown ierr = VecDestroy(&v);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = VecDestroy(&v_vcl);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = VecDestroy(&r_vcl);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = VecDestroy(&d_vcl);CHKERRQ(ierr); 94*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = MatDestroy(&A_vcl);CHKERRQ(ierr); 96*c4762a1bSJed Brown } 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown /* Create a sequential AIJ matrix on rank 0 convert it inplace to a new ViennaCL matrix, and apply it to a ViennaCL vector */ 99*c4762a1bSJed Brown if(!rank) { 100*c4762a1bSJed Brown Mat A; 101*c4762a1bSJed Brown Vec v,r,v_vcl,r_vcl,d_vcl; 102*c4762a1bSJed Brown PetscInt n=17,m=31,nz=5,i,cnt,j; 103*c4762a1bSJed Brown const PetscScalar val = 1.0; 104*c4762a1bSJed Brown PetscReal dnorm; 105*c4762a1bSJed Brown const PetscReal tol = 1e-5; 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,nz,NULL,&A);CHKERRQ(ierr); 108*c4762a1bSJed Brown 109*c4762a1bSJed Brown /* Add nz arbitrary entries per row in arbitrary columns */ 110*c4762a1bSJed Brown for(i=0;i<m;++i) { 111*c4762a1bSJed Brown for(cnt = 0; cnt<nz; ++cnt){ 112*c4762a1bSJed Brown j = (19 * cnt + (7*i + 3)) % n; 113*c4762a1bSJed Brown ierr = MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES);CHKERRQ(ierr); 114*c4762a1bSJed Brown } 115*c4762a1bSJed Brown } 116*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 117*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix");CHKERRQ(ierr); 120*c4762a1bSJed Brown 121*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,n,&v);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,m,&r);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)r,"CPU result vector");CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = VecSet(v,val);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = MatMult(A,v,r);CHKERRQ(ierr); 126*c4762a1bSJed Brown 127*c4762a1bSJed Brown ierr = MatConvert(A,MATSEQAIJVIENNACL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)A,"Converted ViennaCL Matrix");CHKERRQ(ierr); 129*c4762a1bSJed Brown 130*c4762a1bSJed Brown ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl);CHKERRQ(ierr); 131*c4762a1bSJed Brown ierr = VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector");CHKERRQ(ierr); 133*c4762a1bSJed Brown ierr = VecSet(v_vcl,val);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = MatMult(A,v_vcl,r_vcl);CHKERRQ(ierr); 135*c4762a1bSJed Brown 136*c4762a1bSJed Brown ierr = VecDuplicate(r_vcl,&d_vcl);CHKERRQ(ierr); 137*c4762a1bSJed Brown ierr = VecCopy(r_vcl,d_vcl);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = VecAXPY(d_vcl,-1.0,r_vcl); 139*c4762a1bSJed Brown ierr = VecNorm(d_vcl,NORM_INFINITY,&dnorm);CHKERRQ(ierr); 140*c4762a1bSJed Brown if (dnorm > tol) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MPI CPU and MPI ViennaCL Vector results incompatible with inf norm greater than tolerance of %g",tol); 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown ierr = VecDestroy(&v);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 144*c4762a1bSJed Brown ierr = VecDestroy(&v_vcl);CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = VecDestroy(&r_vcl);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = VecDestroy(&d_vcl);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 148*c4762a1bSJed Brown } 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown /* Create a parallel AIJ matrix, convert it to a new ViennaCL matrix, and apply it to a ViennaCL vector */ 151*c4762a1bSJed Brown if (size > 1) { 152*c4762a1bSJed Brown Mat A,A_vcl; 153*c4762a1bSJed Brown Vec v,r,v_vcl,r_vcl,d_vcl; 154*c4762a1bSJed Brown PetscInt N=17,M=31,nz=5,i,cnt,j,rlow,rhigh; 155*c4762a1bSJed Brown const PetscScalar val = 1.0; 156*c4762a1bSJed Brown PetscReal dnorm; 157*c4762a1bSJed Brown const PetscReal tol=1e-5; 158*c4762a1bSJed Brown 159*c4762a1bSJed Brown ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,M,N,nz,NULL,nz,NULL,&A);CHKERRQ(ierr); 160*c4762a1bSJed Brown 161*c4762a1bSJed Brown /* Add nz arbitrary entries per row in arbitrary columns */ 162*c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rlow,&rhigh);CHKERRQ(ierr); 163*c4762a1bSJed Brown for(i=rlow;i<rhigh;++i){ 164*c4762a1bSJed Brown for(cnt = 0; cnt<nz; ++cnt) { 165*c4762a1bSJed Brown j = (19 * cnt + (7*i + 3)) % N; 166*c4762a1bSJed Brown ierr = MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES); 167*c4762a1bSJed Brown } 168*c4762a1bSJed Brown } 169*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 171*c4762a1bSJed Brown 172*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)A,"MPI CPU Matrix");CHKERRQ(ierr); 173*c4762a1bSJed Brown 174*c4762a1bSJed Brown ierr = MatCreateVecs(A,&v,&r);CHKERRQ(ierr); 175*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)r,"MPI CPU result vector");CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = VecSet(v,val);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = MatMult(A,v,r);CHKERRQ(ierr); 178*c4762a1bSJed Brown 179*c4762a1bSJed Brown ierr = MatConvert(A,MATMPIAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl);CHKERRQ(ierr); 180*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)A_vcl,"New MPI ViennaCL Matrix");CHKERRQ(ierr); 181*c4762a1bSJed Brown 182*c4762a1bSJed Brown ierr = MatCreateVecs(A_vcl,&v_vcl,&r_vcl);CHKERRQ(ierr); 183*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector");CHKERRQ(ierr); 184*c4762a1bSJed Brown ierr = VecSet(v_vcl,val);CHKERRQ(ierr); 185*c4762a1bSJed Brown ierr = MatMult(A_vcl,v_vcl,r_vcl);CHKERRQ(ierr); 186*c4762a1bSJed Brown 187*c4762a1bSJed Brown ierr = VecDuplicate(r_vcl,&d_vcl);CHKERRQ(ierr); 188*c4762a1bSJed Brown ierr = VecCopy(r_vcl,d_vcl);CHKERRQ(ierr); 189*c4762a1bSJed Brown ierr = VecAXPY(d_vcl,-1.0,r_vcl); 190*c4762a1bSJed Brown ierr = VecNorm(d_vcl,NORM_INFINITY,&dnorm);CHKERRQ(ierr); 191*c4762a1bSJed Brown if (dnorm > tol) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MPI CPU and MPI ViennaCL Vector results incompatible with inf norm greater than tolerance of %g",tol); 192*c4762a1bSJed Brown 193*c4762a1bSJed Brown ierr = VecDestroy(&v);CHKERRQ(ierr); 194*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 195*c4762a1bSJed Brown ierr = VecDestroy(&v_vcl);CHKERRQ(ierr); 196*c4762a1bSJed Brown ierr = VecDestroy(&r_vcl);CHKERRQ(ierr); 197*c4762a1bSJed Brown ierr = VecDestroy(&d_vcl);CHKERRQ(ierr); 198*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 199*c4762a1bSJed Brown ierr = MatDestroy(&A_vcl);CHKERRQ(ierr); 200*c4762a1bSJed Brown } 201*c4762a1bSJed Brown 202*c4762a1bSJed Brown /* Create a parallel AIJ matrix, convert it in-place to a ViennaCL matrix, and apply it to a ViennaCL vector */ 203*c4762a1bSJed Brown if (size > 1) { 204*c4762a1bSJed Brown Mat A; 205*c4762a1bSJed Brown Vec v,r,v_vcl,r_vcl,d_vcl; 206*c4762a1bSJed Brown PetscInt N=17,M=31,nz=5,i,cnt,j,rlow,rhigh; 207*c4762a1bSJed Brown const PetscScalar val = 1.0; 208*c4762a1bSJed Brown PetscReal dnorm; 209*c4762a1bSJed Brown const PetscReal tol=1e-5; 210*c4762a1bSJed Brown 211*c4762a1bSJed Brown ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,M,N,nz,NULL,nz,NULL,&A);CHKERRQ(ierr); 212*c4762a1bSJed Brown 213*c4762a1bSJed Brown /* Add nz arbitrary entries per row in arbitrary columns */ 214*c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rlow,&rhigh);CHKERRQ(ierr); 215*c4762a1bSJed Brown for(i=rlow;i<rhigh;++i){ 216*c4762a1bSJed Brown for(cnt = 0; cnt<nz; ++cnt) { 217*c4762a1bSJed Brown j = (19 * cnt + (7*i + 3)) % N; 218*c4762a1bSJed Brown ierr = MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES);CHKERRQ(ierr); 219*c4762a1bSJed Brown } 220*c4762a1bSJed Brown } 221*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 222*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 223*c4762a1bSJed Brown 224*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)A,"MPI CPU Matrix");CHKERRQ(ierr); 225*c4762a1bSJed Brown 226*c4762a1bSJed Brown ierr = MatCreateVecs(A,&v,&r);CHKERRQ(ierr); 227*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)r,"MPI CPU result vector");CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = VecSet(v,val);CHKERRQ(ierr); 229*c4762a1bSJed Brown ierr = MatMult(A,v,r);CHKERRQ(ierr); 230*c4762a1bSJed Brown 231*c4762a1bSJed Brown ierr = MatConvert(A,MATMPIAIJVIENNACL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 232*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)A,"Converted MPI ViennaCL Matrix");CHKERRQ(ierr); 233*c4762a1bSJed Brown 234*c4762a1bSJed Brown ierr = MatCreateVecs(A,&v_vcl,&r_vcl);CHKERRQ(ierr); 235*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector");CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = VecSet(v_vcl,val);CHKERRQ(ierr); 237*c4762a1bSJed Brown ierr = MatMult(A,v_vcl,r_vcl);CHKERRQ(ierr); 238*c4762a1bSJed Brown 239*c4762a1bSJed Brown ierr = VecDuplicate(r_vcl,&d_vcl);CHKERRQ(ierr); 240*c4762a1bSJed Brown ierr = VecCopy(r_vcl,d_vcl);CHKERRQ(ierr); 241*c4762a1bSJed Brown ierr = VecAXPY(d_vcl,-1.0,r_vcl);CHKERRQ(ierr); 242*c4762a1bSJed Brown ierr = VecNorm(d_vcl,NORM_INFINITY,&dnorm);CHKERRQ(ierr); 243*c4762a1bSJed Brown if (dnorm > tol) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MPI CPU and MPI ViennaCL Vector results incompatible with inf norm greater than tolerance of %g",tol); 244*c4762a1bSJed Brown 245*c4762a1bSJed Brown ierr = VecDestroy(&v);CHKERRQ(ierr); 246*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 247*c4762a1bSJed Brown ierr = VecDestroy(&v_vcl);CHKERRQ(ierr); 248*c4762a1bSJed Brown ierr = VecDestroy(&r_vcl);CHKERRQ(ierr); 249*c4762a1bSJed Brown ierr = VecDestroy(&d_vcl);CHKERRQ(ierr); 250*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 251*c4762a1bSJed Brown } 252*c4762a1bSJed Brown 253*c4762a1bSJed Brown ierr = PetscFinalize(); 254*c4762a1bSJed Brown return ierr; 255*c4762a1bSJed Brown 256*c4762a1bSJed Brown } 257*c4762a1bSJed Brown 258*c4762a1bSJed Brown /*TEST 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown build: 261*c4762a1bSJed Brown requires: viennacl 262*c4762a1bSJed Brown 263*c4762a1bSJed Brown test: 264*c4762a1bSJed Brown output_file: output/ex204.out 265*c4762a1bSJed Brown 266*c4762a1bSJed Brown test: 267*c4762a1bSJed Brown suffix: 2 268*c4762a1bSJed Brown nsize: 2 269*c4762a1bSJed Brown output_file: output/ex204.out 270*c4762a1bSJed Brown 271*c4762a1bSJed Brown TEST*/ 272