xref: /petsc/src/mat/tests/ex38.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Test interface of Elemental. \n\n";
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            C,Caij;
9*c4762a1bSJed Brown   PetscInt       i,j,m = 5,n,nrows,ncols;
10*c4762a1bSJed Brown   const PetscInt *rows,*cols;
11*c4762a1bSJed Brown   IS             isrows,iscols;
12*c4762a1bSJed Brown   PetscErrorCode ierr;
13*c4762a1bSJed Brown   PetscBool      flg,Test_MatMatMult=PETSC_FALSE,mats_view=PETSC_FALSE;
14*c4762a1bSJed Brown   PetscScalar    *v;
15*c4762a1bSJed Brown   PetscMPIInt    rank,size;
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
18*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
19*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
20*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view);CHKERRQ(ierr);
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown   /* Get local block or element size*/
23*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
24*c4762a1bSJed Brown   n    = m;
25*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = MatSetSizes(C,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = MatSetType(C,MATELEMENTAL);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-row_oriented",&flg);CHKERRQ(ierr);
34*c4762a1bSJed Brown   if (flg) {ierr = MatSetOption(C,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);}
35*c4762a1bSJed Brown   ierr = MatGetOwnershipIS(C,&isrows,&iscols);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-Cexp_view_ownership",&flg);CHKERRQ(ierr);
37*c4762a1bSJed Brown   if (flg) { /* View ownership of explicit C */
38*c4762a1bSJed Brown     IS tmp;
39*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Ownership of explicit C:\n");CHKERRQ(ierr);
40*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Row index set:\n");CHKERRQ(ierr);
41*c4762a1bSJed Brown     ierr = ISOnComm(isrows,PETSC_COMM_WORLD,PETSC_USE_POINTER,&tmp);CHKERRQ(ierr);
42*c4762a1bSJed Brown     ierr = ISView(tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
43*c4762a1bSJed Brown     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
44*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Column index set:\n");CHKERRQ(ierr);
45*c4762a1bSJed Brown     ierr = ISOnComm(iscols,PETSC_COMM_WORLD,PETSC_USE_POINTER,&tmp);CHKERRQ(ierr);
46*c4762a1bSJed Brown     ierr = ISView(tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
47*c4762a1bSJed Brown     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
48*c4762a1bSJed Brown   }
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown   /* Set local matrix entries */
51*c4762a1bSJed Brown   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
52*c4762a1bSJed Brown   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
53*c4762a1bSJed Brown   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
54*c4762a1bSJed Brown   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr);
56*c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
57*c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
58*c4762a1bSJed Brown       /*v[i*ncols+j] = (PetscReal)(rank);*/
59*c4762a1bSJed Brown       v[i*ncols+j] = (PetscReal)(rank*10000+100*rows[i]+cols[j]);
60*c4762a1bSJed Brown     }
61*c4762a1bSJed Brown   }
62*c4762a1bSJed Brown   ierr = MatSetValues(C,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown   /* Test MatView() */
69*c4762a1bSJed Brown   if (mats_view) {
70*c4762a1bSJed Brown     ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
71*c4762a1bSJed Brown   }
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown   /* Test MatMissingDiagonal() */
74*c4762a1bSJed Brown   ierr = MatMissingDiagonal(C,&flg,NULL);CHKERRQ(ierr);
75*c4762a1bSJed Brown   if (flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"MatMissingDiagonal() did not return false!\n");
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown   /* Set unowned matrix entries - add subdiagonals and diagonals from proc[0] */
78*c4762a1bSJed Brown   if (!rank) {
79*c4762a1bSJed Brown     PetscInt M,N,cols[2];
80*c4762a1bSJed Brown     ierr = MatGetSize(C,&M,&N);CHKERRQ(ierr);
81*c4762a1bSJed Brown     for (i=0; i<M; i++) {
82*c4762a1bSJed Brown       cols[0] = i;   v[0] = i + 0.5;
83*c4762a1bSJed Brown       cols[1] = i-1; v[1] = 0.5;
84*c4762a1bSJed Brown       if (i) {
85*c4762a1bSJed Brown         ierr = MatSetValues(C,1,&i,2,cols,v,ADD_VALUES);CHKERRQ(ierr);
86*c4762a1bSJed Brown       } else {
87*c4762a1bSJed Brown         ierr = MatSetValues(C,1,&i,1,&i,v,ADD_VALUES);CHKERRQ(ierr);
88*c4762a1bSJed Brown       }
89*c4762a1bSJed Brown     }
90*c4762a1bSJed Brown   }
91*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
93*c4762a1bSJed Brown 
94*c4762a1bSJed Brown   /* Test MatMult() */
95*c4762a1bSJed Brown   ierr = MatComputeOperator(C,MATAIJ,&Caij);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = MatMultEqual(C,Caij,5,&flg);CHKERRQ(ierr);
97*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultEqual() fails");
98*c4762a1bSJed Brown   ierr = MatMultTransposeEqual(C,Caij,5,&flg);CHKERRQ(ierr);
99*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultTransposeEqual() fails");
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown   /* Test MatMultAdd() and MatMultTransposeAddEqual() */
102*c4762a1bSJed Brown   ierr = MatMultAddEqual(C,Caij,5,&flg);CHKERRQ(ierr);
103*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultAddEqual() fails");
104*c4762a1bSJed Brown   ierr = MatMultTransposeAddEqual(C,Caij,5,&flg);CHKERRQ(ierr);
105*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultTransposeAddEqual() fails");
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown   /* Test MatMatMult() */
108*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-test_matmatmult",&Test_MatMatMult);CHKERRQ(ierr);
109*c4762a1bSJed Brown   if (Test_MatMatMult) {
110*c4762a1bSJed Brown     Mat CCelem,CCaij;
111*c4762a1bSJed Brown     ierr = MatMatMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CCelem);CHKERRQ(ierr);
112*c4762a1bSJed Brown     ierr = MatMatMult(Caij,Caij,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CCaij);CHKERRQ(ierr);
113*c4762a1bSJed Brown     ierr = MatMultEqual(CCelem,CCaij,5,&flg);CHKERRQ(ierr);
114*c4762a1bSJed Brown     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"CCelem != CCaij. MatMatMult() fails");
115*c4762a1bSJed Brown     ierr = MatDestroy(&CCaij);CHKERRQ(ierr);
116*c4762a1bSJed Brown     ierr = MatDestroy(&CCelem);CHKERRQ(ierr);
117*c4762a1bSJed Brown   }
118*c4762a1bSJed Brown 
119*c4762a1bSJed Brown   ierr = PetscFree(v);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = MatDestroy(&Caij);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = PetscFinalize();
125*c4762a1bSJed Brown   return ierr;
126*c4762a1bSJed Brown }
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown 
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown /*TEST
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown    test:
134*c4762a1bSJed Brown       nsize: 2
135*c4762a1bSJed Brown       requires: elemental
136*c4762a1bSJed Brown       args: -mat_type elemental -m 2 -n 3
137*c4762a1bSJed Brown 
138*c4762a1bSJed Brown    test:
139*c4762a1bSJed Brown       suffix: 2
140*c4762a1bSJed Brown       nsize: 6
141*c4762a1bSJed Brown       requires: elemental
142*c4762a1bSJed Brown       args: -mat_type elemental -m 2 -n 2
143*c4762a1bSJed Brown 
144*c4762a1bSJed Brown    test:
145*c4762a1bSJed Brown       suffix: 3
146*c4762a1bSJed Brown       nsize: 6
147*c4762a1bSJed Brown       requires: elemental
148*c4762a1bSJed Brown       args: -mat_type elemental -m 2 -n 2 -test_matmatmult
149*c4762a1bSJed Brown       output_file: output/ex38_2.out
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown TEST*/
152