xref: /petsc/src/mat/tests/ex26.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests MatGetRowIJ for SeqAIJ, SeqBAIJ and SeqSBAIJ\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown PetscErrorCode DumpCSR(Mat A,PetscInt shift,PetscBool symmetric,PetscBool compressed)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   MatType        type;
9*c4762a1bSJed Brown   PetscInt       i,j,nr,bs = 1;
10*c4762a1bSJed Brown   const PetscInt *ia,*ja;
11*c4762a1bSJed Brown   PetscBool      done,isseqbaij,isseqsbaij;
12*c4762a1bSJed Brown   PetscErrorCode ierr;
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown   PetscFunctionBeginUser;
15*c4762a1bSJed Brown   ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQBAIJ,&isseqbaij);CHKERRQ(ierr);
16*c4762a1bSJed Brown   ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
17*c4762a1bSJed Brown   if (isseqbaij || isseqsbaij) {
18*c4762a1bSJed Brown     ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
19*c4762a1bSJed Brown   }
20*c4762a1bSJed Brown   ierr = MatGetType(A,&type);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = MatGetRowIJ(A,shift,symmetric,compressed,&nr,&ia,&ja,&done);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"===========================================================\n");CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"CSR for %s: shift %D symmetric %D compressed %D\n",type,shift,(PetscInt)symmetric,(PetscInt)compressed);CHKERRQ(ierr);
24*c4762a1bSJed Brown   for (i=0;i<nr;i++) {
25*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"%D:",i+shift);CHKERRQ(ierr);
26*c4762a1bSJed Brown     for (j=ia[i];j<ia[i+1];j++) {
27*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF," %D",ja[j-shift]);CHKERRQ(ierr);
28*c4762a1bSJed Brown     }
29*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"\n");CHKERRQ(ierr);
30*c4762a1bSJed Brown   }
31*c4762a1bSJed Brown   ierr = MatRestoreRowIJ(A,shift,symmetric,compressed,&nr,&ia,&ja,&done);CHKERRQ(ierr);
32*c4762a1bSJed Brown   PetscFunctionReturn(0);
33*c4762a1bSJed Brown }
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown int main(int argc,char **args)
36*c4762a1bSJed Brown {
37*c4762a1bSJed Brown   Mat            A,B,C;
38*c4762a1bSJed Brown   PetscInt       i,j,k,m = 3,n = 3,bs = 1;
39*c4762a1bSJed Brown   PetscErrorCode ierr;
40*c4762a1bSJed Brown   PetscMPIInt    size;
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
43*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
44*c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
45*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr);
48*c4762a1bSJed Brown   /* adjust sizes by block size */
49*c4762a1bSJed Brown   if (m%bs) m += bs-m%bs;
50*c4762a1bSJed Brown   if (n%bs) n += bs-n%bs;
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr);
53*c4762a1bSJed Brown   ierr = MatSetSizes(A,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
54*c4762a1bSJed Brown   ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = MatSetSizes(B,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = MatSetBlockSize(B,bs);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = MatSetType(B,MATSEQBAIJ);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = MatSetUp(B);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&C);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = MatSetSizes(C,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = MatSetBlockSize(C,bs);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown   for (i=0; i<m; i++) {
70*c4762a1bSJed Brown     for (j=0; j<n; j++) {
71*c4762a1bSJed Brown 
72*c4762a1bSJed Brown       PetscScalar v = -1.0;
73*c4762a1bSJed Brown       PetscInt    Ii = j + n*i,J;
74*c4762a1bSJed Brown       J = Ii - n;
75*c4762a1bSJed Brown       if (J>=0)  {
76*c4762a1bSJed Brown         ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
77*c4762a1bSJed Brown         ierr = MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
78*c4762a1bSJed Brown         ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
79*c4762a1bSJed Brown       }
80*c4762a1bSJed Brown       J = Ii + n;
81*c4762a1bSJed Brown       if (J<m*n) {
82*c4762a1bSJed Brown         ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
83*c4762a1bSJed Brown         ierr = MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
84*c4762a1bSJed Brown         ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
85*c4762a1bSJed Brown       }
86*c4762a1bSJed Brown       J = Ii - 1;
87*c4762a1bSJed Brown       if (J>=0)  {
88*c4762a1bSJed Brown         ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
89*c4762a1bSJed Brown         ierr = MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
90*c4762a1bSJed Brown         ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
91*c4762a1bSJed Brown       }
92*c4762a1bSJed Brown       J = Ii + 1;
93*c4762a1bSJed Brown       if (J<m*n) {
94*c4762a1bSJed Brown         ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
95*c4762a1bSJed Brown         ierr = MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
96*c4762a1bSJed Brown         ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
97*c4762a1bSJed Brown       }
98*c4762a1bSJed Brown       v = 4.0;
99*c4762a1bSJed Brown       ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
100*c4762a1bSJed Brown       ierr = MatSetValues(B,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
101*c4762a1bSJed Brown       ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
102*c4762a1bSJed Brown     }
103*c4762a1bSJed Brown   }
104*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown   /* test MatGetRowIJ for the three Mat types */
112*c4762a1bSJed Brown   ierr = MatView(A,NULL);CHKERRQ(ierr);
113*c4762a1bSJed Brown   ierr = MatView(B,NULL);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = MatView(C,NULL);CHKERRQ(ierr);
115*c4762a1bSJed Brown   for (i=0;i<2;i++) {
116*c4762a1bSJed Brown     PetscInt shift = i;
117*c4762a1bSJed Brown     for (j=0;j<2;j++) {
118*c4762a1bSJed Brown       PetscBool symmetric = ((j>0) ? PETSC_FALSE : PETSC_TRUE);
119*c4762a1bSJed Brown       for (k=0;k<2;k++) {
120*c4762a1bSJed Brown         PetscBool compressed = ((k>0) ? PETSC_FALSE : PETSC_TRUE);
121*c4762a1bSJed Brown         ierr = DumpCSR(A,shift,symmetric,compressed);CHKERRQ(ierr);
122*c4762a1bSJed Brown         ierr = DumpCSR(B,shift,symmetric,compressed);CHKERRQ(ierr);
123*c4762a1bSJed Brown         ierr = DumpCSR(C,shift,symmetric,compressed);CHKERRQ(ierr);
124*c4762a1bSJed Brown       }
125*c4762a1bSJed Brown     }
126*c4762a1bSJed Brown   }
127*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
129*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
130*c4762a1bSJed Brown   ierr = PetscFinalize();
131*c4762a1bSJed Brown   return ierr;
132*c4762a1bSJed Brown }
133*c4762a1bSJed Brown 
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown /*TEST
136*c4762a1bSJed Brown 
137*c4762a1bSJed Brown    test:
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown    test:
140*c4762a1bSJed Brown       suffix: 2
141*c4762a1bSJed Brown       args: -bs 2
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown TEST*/
144