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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 10c4762a1bSJed Brown 115f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 125f80ce2aSJacob Faibussowitsch CHKERRMPI(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 215f80ce2aSJacob Faibussowitsch CHKERRQ(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; 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A_vcl,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 28c4762a1bSJed Brown } 29c4762a1bSJed Brown } 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A_vcl,MAT_FINAL_ASSEMBLY)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A_vcl,MAT_FINAL_ASSEMBLY)); 32c4762a1bSJed Brown 335f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_vcl,val)); 36c4762a1bSJed Brown 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A_vcl,v_vcl,r_vcl)); 38c4762a1bSJed Brown 395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_vcl)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r_vcl)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(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 535f80ce2aSJacob Faibussowitsch CHKERRQ(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; 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,i,j,(PetscScalar) (0.3 * i + j),INSERT_VALUES)); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown } 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 64c4762a1bSJed Brown 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix")); 66c4762a1bSJed Brown 675f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&v)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&r)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r,"CPU result vector")); 705f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,val)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v,r)); 72c4762a1bSJed Brown 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSEQAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A_vcl,"New ViennaCL Matrix")); 75c4762a1bSJed Brown 765f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 795f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_vcl,val)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A_vcl,v_vcl,r_vcl)); 81c4762a1bSJed Brown 825f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(r_vcl,&d_vcl)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(r_vcl,d_vcl)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(d_vcl,-1.0,r_vcl)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(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 885f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_vcl)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r_vcl)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d_vcl)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(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 1065f80ce2aSJacob Faibussowitsch CHKERRQ(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; 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown } 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 117c4762a1bSJed Brown 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"Sequential CPU Matrix")); 119c4762a1bSJed Brown 1205f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&v)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&r)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r,"CPU result vector")); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,val)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v,r)); 125c4762a1bSJed Brown 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSEQAIJVIENNACL,MAT_INPLACE_MATRIX,&A)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"Converted ViennaCL Matrix")); 128c4762a1bSJed Brown 1295f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,n,&v_vcl)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&r_vcl)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_vcl,val)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v_vcl,r_vcl)); 134c4762a1bSJed Brown 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(r_vcl,&d_vcl)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(r_vcl,d_vcl)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(d_vcl,-1.0,r_vcl)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(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 1415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_vcl)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r_vcl)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d_vcl)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(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 1585f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 1615f80ce2aSJacob Faibussowitsch CHKERRQ(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; 1655f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown } 1685f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 170c4762a1bSJed Brown 1715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"MPI CPU Matrix")); 172c4762a1bSJed Brown 1735f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&v,&r)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r,"MPI CPU result vector")); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,val)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v,r)); 177c4762a1bSJed Brown 1785f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATMPIAIJVIENNACL,MAT_INITIAL_MATRIX,&A_vcl)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A_vcl,"New MPI ViennaCL Matrix")); 180c4762a1bSJed Brown 1815f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A_vcl,&v_vcl,&r_vcl)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_vcl,val)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A_vcl,v_vcl,r_vcl)); 185c4762a1bSJed Brown 1865f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(r_vcl,&d_vcl)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(r_vcl,d_vcl)); 1885f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(d_vcl,-1.0,r_vcl)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(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 1925f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_vcl)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r_vcl)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d_vcl)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(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 2105f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 2135f80ce2aSJacob Faibussowitsch CHKERRQ(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; 2175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,i,j,(PetscScalar)(0.3 * i + j),INSERT_VALUES)); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown } 2205f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 222c4762a1bSJed Brown 2235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"MPI CPU Matrix")); 224c4762a1bSJed Brown 2255f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&v,&r)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r,"MPI CPU result vector")); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,val)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v,r)); 229c4762a1bSJed Brown 2305f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATMPIAIJVIENNACL,MAT_INPLACE_MATRIX,&A)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"Converted MPI ViennaCL Matrix")); 232c4762a1bSJed Brown 2335f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&v_vcl,&r_vcl)); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)r_vcl,"ViennaCL result vector")); 2355f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_vcl,val)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v_vcl,r_vcl)); 237c4762a1bSJed Brown 2385f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(r_vcl,&d_vcl)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(r_vcl,d_vcl)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(d_vcl,-1.0,r_vcl)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(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 2445f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_vcl)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r_vcl)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d_vcl)); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 253*b122ec5aSJacob 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