1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatGetColumnVector()."; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat C; 9c4762a1bSJed Brown PetscInt i,j,m = 3,n = 2,Ii,J,col = 0; 10c4762a1bSJed Brown PetscMPIInt size,rank; 11c4762a1bSJed Brown PetscErrorCode ierr; 12c4762a1bSJed Brown PetscScalar v; 13c4762a1bSJed Brown Vec yy; 14c4762a1bSJed Brown 15c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 16*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-col",&col,NULL)); 17c4762a1bSJed Brown 18*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 19*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 20c4762a1bSJed Brown n = 2*size; 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN*/ 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(C,5,NULL)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(C,5,NULL,5,NULL)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 29c4762a1bSJed Brown 30c4762a1bSJed Brown for (i=0; i<m; i++) { 31c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 32c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 33*5f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 34*5f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 35*5f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 36*5f80ce2aSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 37*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown } 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 43c4762a1bSJed Brown 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(C,NULL,&yy)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(yy)); 46c4762a1bSJed Brown 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetColumnVector(C,yy,col)); 48c4762a1bSJed Brown 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(yy,PETSC_VIEWER_STDOUT_WORLD)); 50c4762a1bSJed Brown 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&yy)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 53c4762a1bSJed Brown ierr = PetscFinalize(); 54c4762a1bSJed Brown return ierr; 55c4762a1bSJed Brown } 56c4762a1bSJed Brown 57c4762a1bSJed Brown /*TEST 58c4762a1bSJed Brown 59c4762a1bSJed Brown test: 60c4762a1bSJed Brown nsize: 3 61c4762a1bSJed Brown args: -col 7 62c4762a1bSJed Brown 63c4762a1bSJed Brown test: 64c4762a1bSJed Brown suffix: dense 65c4762a1bSJed Brown nsize: 3 66a5225ed3SStefano Zampini args: -col 7 -mat_type dense -vec_type {{mpi standard}} 67a5225ed3SStefano Zampini filter: grep -v type 68a5225ed3SStefano Zampini 69a5225ed3SStefano Zampini test: 70a5225ed3SStefano Zampini requires: cuda 71a5225ed3SStefano Zampini suffix: dense_cuda 72a5225ed3SStefano Zampini nsize: 3 73a5225ed3SStefano Zampini output_file: output/ex60_dense.out 74a5225ed3SStefano Zampini args: -col 7 -mat_type {{mpidense mpidensecuda}} -vec_type {{mpi standard cuda mpicuda}} 75a5225ed3SStefano Zampini filter: grep -v type 76c4762a1bSJed Brown 77c4762a1bSJed Brown TEST*/ 78