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