xref: /petsc/src/mat/tests/ex60.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscScalar    v;
12c4762a1bSJed Brown   Vec            yy;
13c4762a1bSJed Brown 
14*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-col",&col,NULL));
16c4762a1bSJed Brown 
175f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
185f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
19c4762a1bSJed Brown   n    = 2*size;
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(C,5,NULL));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(C,5,NULL,5,NULL));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   for (i=0; i<m; i++) {
30c4762a1bSJed Brown     for (j=2*rank; j<2*rank+2; j++) {
31c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
325f80ce2aSJacob Faibussowitsch       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
335f80ce2aSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
345f80ce2aSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
355f80ce2aSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
365f80ce2aSJacob Faibussowitsch       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
37c4762a1bSJed Brown     }
38c4762a1bSJed Brown   }
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
42c4762a1bSJed Brown 
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(C,NULL,&yy));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(yy));
45c4762a1bSJed Brown 
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetColumnVector(C,yy,col));
47c4762a1bSJed Brown 
485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(yy,PETSC_VIEWER_STDOUT_WORLD));
49c4762a1bSJed Brown 
505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&yy));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
52*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
53*b122ec5aSJacob Faibussowitsch   return 0;
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /*TEST
57c4762a1bSJed Brown 
58c4762a1bSJed Brown    test:
59c4762a1bSJed Brown       nsize: 3
60c4762a1bSJed Brown       args: -col 7
61c4762a1bSJed Brown 
62c4762a1bSJed Brown    test:
63c4762a1bSJed Brown       suffix: dense
64c4762a1bSJed Brown       nsize: 3
65a5225ed3SStefano Zampini       args: -col 7 -mat_type dense -vec_type {{mpi standard}}
66a5225ed3SStefano Zampini       filter: grep -v type
67a5225ed3SStefano Zampini 
68a5225ed3SStefano Zampini    test:
69a5225ed3SStefano Zampini       requires: cuda
70a5225ed3SStefano Zampini       suffix: dense_cuda
71a5225ed3SStefano Zampini       nsize: 3
72a5225ed3SStefano Zampini       output_file: output/ex60_dense.out
73a5225ed3SStefano Zampini       args: -col 7 -mat_type {{mpidense mpidensecuda}} -vec_type {{mpi standard cuda mpicuda}}
74a5225ed3SStefano Zampini       filter: grep -v type
75c4762a1bSJed Brown 
76c4762a1bSJed Brown TEST*/
77