xref: /petsc/src/mat/tests/ex52.c (revision ffc4695bcb29f4b022f59a5fd6bc99fc280ff6d8)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPIBAIJ format.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscmat.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            A;
10c4762a1bSJed Brown   PetscInt       m=2,bs=1,M,row,col,start,end,i,j,k;
11c4762a1bSJed Brown   PetscErrorCode ierr;
12c4762a1bSJed Brown   PetscMPIInt    rank,size;
13c4762a1bSJed Brown   PetscScalar    data=100;
14c4762a1bSJed Brown   PetscBool      flg;
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
17*ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
18*ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   /* Test MatSetValues() and MatGetValues() */
21c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL);CHKERRQ(ierr);
22c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL);CHKERRQ(ierr);
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   M    = m*bs*size;
25c4762a1bSJed Brown   ierr = MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,M,M,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A);CHKERRQ(ierr);
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&start,&end);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-column_oriented",&flg);CHKERRQ(ierr);
29c4762a1bSJed Brown   if (flg) {
30c4762a1bSJed Brown     ierr = MatSetOption(A,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
31c4762a1bSJed Brown   }
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /* inproc assembly */
34c4762a1bSJed Brown   for (row=start; row<end; row++) {
35c4762a1bSJed Brown     for (col=start; col<end; col++,data+=1) {
36c4762a1bSJed Brown       ierr = MatSetValues(A,1,&row,1,&col,&data,INSERT_VALUES);CHKERRQ(ierr);
37c4762a1bSJed Brown     }
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* offproc assembly */
43c4762a1bSJed Brown   data = 5.0;
44c4762a1bSJed Brown   row  = (M+start-1)%M;
45c4762a1bSJed Brown   for (col=0; col<M; col++) {
46c4762a1bSJed Brown     ierr = MatSetValues(A,1,&row,1,&col,&data,ADD_VALUES);CHKERRQ(ierr);
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Test MatSetValuesBlocked() */
52c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-test_setvaluesblocked",&flg);CHKERRQ(ierr);
53c4762a1bSJed Brown   if (flg) {
54c4762a1bSJed Brown     PetscScalar *bval;
55c4762a1bSJed Brown     row /= bs;
56c4762a1bSJed Brown     col  = start/bs;
57c4762a1bSJed Brown     ierr = PetscMalloc1(bs*bs,&bval);CHKERRQ(ierr);
58c4762a1bSJed Brown     k = 1;
59c4762a1bSJed Brown     /* row oriented - defalt */
60c4762a1bSJed Brown     for (i=0; i<bs; i++) {
61c4762a1bSJed Brown       for (j=0; j<bs; j++) {
62c4762a1bSJed Brown         bval[i*bs+j] = (PetscScalar)k; k++;
63c4762a1bSJed Brown       }
64c4762a1bSJed Brown     }
65c4762a1bSJed Brown     ierr = MatSetValuesBlocked(A,1,&row,1,&col,bval,INSERT_VALUES);CHKERRQ(ierr);
66c4762a1bSJed Brown     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67c4762a1bSJed Brown     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68c4762a1bSJed Brown     ierr = PetscFree(bval);CHKERRQ(ierr);
69c4762a1bSJed Brown   }
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
72c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
73c4762a1bSJed Brown   ierr = PetscFinalize();
74c4762a1bSJed Brown   return ierr;
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown 
78c4762a1bSJed Brown 
79c4762a1bSJed Brown /*TEST
80c4762a1bSJed Brown 
81c4762a1bSJed Brown    test:
82c4762a1bSJed Brown       suffix: 1
83c4762a1bSJed Brown       nsize: 3
84c4762a1bSJed Brown       args: -mat_block_size 2 -test_setvaluesblocked
85c4762a1bSJed Brown 
86c4762a1bSJed Brown    test:
87c4762a1bSJed Brown       suffix: 2
88c4762a1bSJed Brown       nsize: 3
89c4762a1bSJed Brown       args: -mat_block_size 2 -test_setvaluesblocked -column_oriented
90c4762a1bSJed Brown 
91c4762a1bSJed Brown    test:
92c4762a1bSJed Brown       suffix: 3
93c4762a1bSJed Brown       nsize: 3
94c4762a1bSJed Brown       args: -mat_block_size 1 -test_setvaluesblocked
95c4762a1bSJed Brown 
96c4762a1bSJed Brown    test:
97c4762a1bSJed Brown       suffix: 4
98c4762a1bSJed Brown       nsize: 3
99c4762a1bSJed Brown       args: -mat_block_size 1 -test_setvaluesblocked -column_oriented
100c4762a1bSJed Brown 
101c4762a1bSJed Brown TEST*/
102c4762a1bSJed Brown 
103