xref: /petsc/src/mat/tests/ex56.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Test the use of MatSetValuesBlocked(), MatZeroRows() for rectangular MatBAIJ matrix, test MatSetValuesBlocked() for MatSBAIJ matrix (-test_mat_sbaij).";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **args)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            A;
9*c4762a1bSJed Brown   PetscInt       bs=3,m=4,n=6,i,j,val = 10,row[2],col[3],eval,rstart;
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   PetscMPIInt    size,rank;
12*c4762a1bSJed Brown   PetscScalar    x[6][9],y[3][3],one=1.0;
13*c4762a1bSJed Brown   PetscBool      flg,testsbaij=PETSC_FALSE;
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
16*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
17*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-test_mat_sbaij",&testsbaij);CHKERRQ(ierr);
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown   if (testsbaij) {
22*c4762a1bSJed Brown     ierr = MatCreateSBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,NULL,1,NULL,&A);CHKERRQ(ierr);
23*c4762a1bSJed Brown   } else {
24*c4762a1bSJed Brown     ierr = MatCreateBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,NULL,1,NULL,&A);CHKERRQ(ierr);
25*c4762a1bSJed Brown   }
26*c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
27*c4762a1bSJed Brown   eval = 9;
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-ass_extern",&flg);CHKERRQ(ierr);
30*c4762a1bSJed Brown   if (flg && (size != 1)) rstart = m*((rank+1)%size);
31*c4762a1bSJed Brown   else rstart = m*(rank);
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   row[0] =rstart+0;  row[1] =rstart+2;
34*c4762a1bSJed Brown   col[0] =rstart+0;  col[1] =rstart+1;  col[2] =rstart+3;
35*c4762a1bSJed Brown   for (i=0; i<6; i++) {
36*c4762a1bSJed Brown     for (j =0; j< 9; j++) x[i][j] = (PetscScalar)val++;
37*c4762a1bSJed Brown   }
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown   ierr = MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES);CHKERRQ(ierr);
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown   /*
47*c4762a1bSJed Brown   This option does not work for rectangular matrices
48*c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
49*c4762a1bSJed Brown   */
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown   ierr = MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES);CHKERRQ(ierr);
52*c4762a1bSJed Brown 
53*c4762a1bSJed Brown   /* Do another MatSetValues to test the case when only one local block is specified */
54*c4762a1bSJed Brown   for (i=0; i<3; i++) {
55*c4762a1bSJed Brown     for (j =0; j<3; j++) y[i][j] = (PetscScalar)(10 + i*eval + j);
56*c4762a1bSJed Brown   }
57*c4762a1bSJed Brown   ierr = MatSetValuesBlocked(A,1,row,1,col,&y[0][0],INSERT_VALUES);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60*c4762a1bSJed Brown 
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-zero_rows",&flg);CHKERRQ(ierr);
63*c4762a1bSJed Brown   if (flg) {
64*c4762a1bSJed Brown     col[0] = rstart*bs+0;
65*c4762a1bSJed Brown     col[1] = rstart*bs+1;
66*c4762a1bSJed Brown     col[2] = rstart*bs+2;
67*c4762a1bSJed Brown     ierr   = MatZeroRows(A,3,col,one,0,0);CHKERRQ(ierr);
68*c4762a1bSJed Brown   }
69*c4762a1bSJed Brown 
70*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
71*c4762a1bSJed Brown 
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 /*TEST
79*c4762a1bSJed Brown 
80*c4762a1bSJed Brown    test:
81*c4762a1bSJed Brown       filter: grep -v "MPI processes"
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown    test:
84*c4762a1bSJed Brown       suffix: 4
85*c4762a1bSJed Brown       nsize: 3
86*c4762a1bSJed Brown       args: -ass_extern
87*c4762a1bSJed Brown       filter: grep -v "MPI processes"
88*c4762a1bSJed Brown 
89*c4762a1bSJed Brown    test:
90*c4762a1bSJed Brown       suffix: 5
91*c4762a1bSJed Brown       nsize: 3
92*c4762a1bSJed Brown       args: -ass_extern -zero_rows
93*c4762a1bSJed Brown       filter: grep -v "MPI processes"
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown TEST*/
96