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