xref: /petsc/src/mat/tests/ex300.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Show MatShift BUG happening after copying a matrix with no rows on a process";
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    Contributed by: Eric Chamberland
5c4762a1bSJed Brown */
6c4762a1bSJed Brown #include <petscmat.h>
7c4762a1bSJed Brown 
8c4762a1bSJed Brown /* DEFINE this to turn on/off the bug: */
9c4762a1bSJed Brown #define SET_2nd_PROC_TO_HAVE_NO_LOCAL_LINES
10c4762a1bSJed Brown 
11c4762a1bSJed Brown int main(int argc,char **args)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   Mat               C;
14c4762a1bSJed Brown   PetscInt          i,m = 3;
15c4762a1bSJed Brown   PetscMPIInt       rank,size;
16c4762a1bSJed Brown   PetscErrorCode    ierr;
17c4762a1bSJed Brown   PetscScalar       v;
18c4762a1bSJed Brown   Mat               lMatA;
19c4762a1bSJed Brown   PetscInt          locallines;
20c4762a1bSJed Brown   PetscInt          d_nnz[3] = {0,0,0};
21c4762a1bSJed Brown   PetscInt          o_nnz[3] = {0,0,0};
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
24ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
25ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
26c4762a1bSJed Brown 
27*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(2 != size,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"Relevant with 2 processes only");
28c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
29c4762a1bSJed Brown 
30c4762a1bSJed Brown #ifdef SET_2nd_PROC_TO_HAVE_NO_LOCAL_LINES
31c4762a1bSJed Brown   if (0 == rank) {
32c4762a1bSJed Brown     locallines = m;
33c4762a1bSJed Brown     d_nnz[0] = 1;
34c4762a1bSJed Brown     d_nnz[1] = 1;
35c4762a1bSJed Brown     d_nnz[2] = 1;
36c4762a1bSJed Brown   } else {
37c4762a1bSJed Brown    locallines = 0;
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown #else
40c4762a1bSJed Brown   if (0 == rank) {
41c4762a1bSJed Brown     locallines = m-1;
42c4762a1bSJed Brown     d_nnz[0] = 1;
43c4762a1bSJed Brown     d_nnz[1] = 1;
44c4762a1bSJed Brown   } else {
45c4762a1bSJed Brown     locallines = 1;
46c4762a1bSJed Brown     d_nnz[0] = 1;
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown #endif
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   ierr = MatSetSizes(C,locallines,locallines,m,m);CHKERRQ(ierr);
51c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = MatXAIJSetPreallocation(C,1,d_nnz,o_nnz,NULL,NULL);CHKERRQ(ierr);
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   v = 2;
55c4762a1bSJed Brown   /* Assembly on the diagonal: */
56c4762a1bSJed Brown   for (i=0; i<m; i++) {
57c4762a1bSJed Brown      ierr = MatSetValues(C,1,&i,1,&i,&v,ADD_VALUES);CHKERRQ(ierr);
58c4762a1bSJed Brown   }
59c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
61c4762a1bSJed Brown   ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr = MatSetOption(C, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE);CHKERRQ(ierr);
63c4762a1bSJed Brown   ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
64c4762a1bSJed Brown   ierr = MatConvert(C,MATSAME, MAT_INITIAL_MATRIX, &lMatA);CHKERRQ(ierr);
65c4762a1bSJed Brown   ierr = MatView(lMatA,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   ierr = MatShift(lMatA,-1.0);CHKERRQ(ierr);
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   ierr = MatDestroy(&lMatA);CHKERRQ(ierr);
70c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
71c4762a1bSJed Brown   ierr = PetscFinalize();
72c4762a1bSJed Brown   return ierr;
73c4762a1bSJed Brown }
74c4762a1bSJed Brown 
75c4762a1bSJed Brown /*TEST
76c4762a1bSJed Brown 
77c4762a1bSJed Brown    test:
78c4762a1bSJed Brown       nsize: 2
79c4762a1bSJed Brown 
80c4762a1bSJed Brown TEST*/
81