xref: /petsc/src/mat/tests/ex60.c (revision a5225ed31aa9ee76d5b3e7690f8edc53c82a15a1)
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;
16c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-col",&col,NULL);CHKERRQ(ierr);
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
19c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
20c4762a1bSJed Brown   n    = 2*size;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
23c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
24c4762a1bSJed Brown   ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
26c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(C,5,NULL);CHKERRQ(ierr);
27c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocation(C,5,NULL,5,NULL);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
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;
33c4762a1bSJed Brown       if (i>0)   {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
34c4762a1bSJed Brown       if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
35c4762a1bSJed Brown       if (j>0)   {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
36c4762a1bSJed Brown       if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
37c4762a1bSJed Brown       v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
38c4762a1bSJed Brown     }
39c4762a1bSJed Brown   }
40c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
43c4762a1bSJed Brown 
44*a5225ed3SStefano Zampini   ierr = MatCreateVecs(C,NULL,&yy);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = VecSetFromOptions(yy);CHKERRQ(ierr);
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   ierr = MatGetColumnVector(C,yy,col);CHKERRQ(ierr);
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   ierr = VecView(yy,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   ierr = VecDestroy(&yy);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
53c4762a1bSJed Brown   ierr = PetscFinalize();
54c4762a1bSJed Brown   return ierr;
55c4762a1bSJed Brown }
56c4762a1bSJed Brown 
57c4762a1bSJed Brown 
58c4762a1bSJed Brown /*TEST
59c4762a1bSJed Brown 
60c4762a1bSJed Brown    test:
61c4762a1bSJed Brown       nsize: 3
62c4762a1bSJed Brown       args: -col 7
63c4762a1bSJed Brown 
64c4762a1bSJed Brown    test:
65c4762a1bSJed Brown       suffix: dense
66c4762a1bSJed Brown       nsize: 3
67*a5225ed3SStefano Zampini       args: -col 7 -mat_type dense -vec_type {{mpi standard}}
68*a5225ed3SStefano Zampini       filter: grep -v type
69*a5225ed3SStefano Zampini 
70*a5225ed3SStefano Zampini    test:
71*a5225ed3SStefano Zampini       requires: cuda
72*a5225ed3SStefano Zampini       suffix: dense_cuda
73*a5225ed3SStefano Zampini       nsize: 3
74*a5225ed3SStefano Zampini       output_file: output/ex60_dense.out
75*a5225ed3SStefano Zampini       args: -col 7 -mat_type {{mpidense mpidensecuda}} -vec_type {{mpi standard cuda mpicuda}}
76*a5225ed3SStefano Zampini       filter: grep -v type
77c4762a1bSJed Brown 
78c4762a1bSJed Brown TEST*/
79