xref: /petsc/src/mat/tests/ex60.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry Smith   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-col",&col,NULL));
17c4762a1bSJed Brown 
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
199566063dSJacob Faibussowitsch   PetscCallMPI(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*/
239566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
249566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
259566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
269566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(C,5,NULL));
279566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(C,5,NULL,5,NULL));
289566063dSJacob Faibussowitsch   PetscCall(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;
339566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
349566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
359566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
369566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
379566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
38c4762a1bSJed Brown     }
39c4762a1bSJed Brown   }
409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
429566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
43c4762a1bSJed Brown 
449566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(C,NULL,&yy));
459566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(yy));
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(MatGetColumnVector(C,yy,col));
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCall(VecView(yy,PETSC_VIEWER_STDOUT_WORLD));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&yy));
529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
539566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
54b122ec5aSJacob Faibussowitsch   return 0;
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