xref: /petsc/src/mat/tests/ex246.cxx (revision 70a7d78aacfbd24b2e31399a3d8e056944bb7de3)
1c7a4214aSPierre Jolivet static char help[] = "Tests MATHTOOL with a derived htool::IMatrix<PetscScalar> class\n\n";
2c7a4214aSPierre Jolivet 
3c7a4214aSPierre Jolivet #include <petscmat.h>
4c7a4214aSPierre Jolivet #include <htool/misc/petsc.hpp>
5c7a4214aSPierre Jolivet 
6c7a4214aSPierre Jolivet static PetscScalar GenEntry(PetscInt sdim,PetscInt i,PetscInt j,void *ctx)
7c7a4214aSPierre Jolivet {
8c7a4214aSPierre Jolivet   PetscInt  d;
9c7a4214aSPierre Jolivet   PetscReal diff = 0.0,*coords = (PetscReal*)(ctx);
10c7a4214aSPierre Jolivet 
11c7a4214aSPierre Jolivet   PetscFunctionBeginUser;
12c7a4214aSPierre Jolivet   for (d = 0; d < sdim; d++) { diff += (coords[i*sdim+d] - coords[j*sdim+d]) * (coords[i*sdim+d] - coords[j*sdim+d]); }
13c7a4214aSPierre Jolivet   PetscFunctionReturn(1.0/(1.0e-2 + PetscSqrtReal(diff)));
14c7a4214aSPierre Jolivet }
15c7a4214aSPierre Jolivet 
16c7a4214aSPierre Jolivet class MyIMatrix : public htool::IMatrix<PetscScalar> {
17c7a4214aSPierre Jolivet   private:
18c7a4214aSPierre Jolivet   PetscReal *coords;
19c7a4214aSPierre Jolivet   PetscInt  sdim;
20c7a4214aSPierre Jolivet   public:
21c7a4214aSPierre Jolivet   MyIMatrix(PetscInt M,PetscInt N,PetscInt spacedim,PetscReal* gcoords) : htool::IMatrix<PetscScalar>(M,N),coords(gcoords),sdim(spacedim) { }
22c7a4214aSPierre Jolivet 
23*70a7d78aSStefano Zampini   PetscScalar get_coef(const PetscInt &i, const PetscInt &j) const override
24*70a7d78aSStefano Zampini   {
25c7a4214aSPierre Jolivet     PetscReal diff = 0.0;
26c7a4214aSPierre Jolivet 
27c7a4214aSPierre Jolivet     PetscFunctionBeginUser;
28c7a4214aSPierre Jolivet     for (PetscInt d = 0; d < sdim; d++) { diff += (coords[i*sdim+d] - coords[j*sdim+d]) * (coords[i*sdim+d] - coords[j*sdim+d]); }
29c7a4214aSPierre Jolivet     PetscFunctionReturn(1.0/(1.0e-2 + PetscSqrtReal(diff)));
30c7a4214aSPierre Jolivet   }
31c7a4214aSPierre Jolivet 
32*70a7d78aSStefano Zampini   void copy_submatrix(PetscInt M, PetscInt N, const PetscInt *const rows, const PetscInt *const cols, PetscScalar *ptr) const override
33*70a7d78aSStefano Zampini   {
34c7a4214aSPierre Jolivet     PetscFunctionBeginUser;
35c7a4214aSPierre Jolivet     for (PetscInt j = 0; j < M; j++) /* could be optimized by the user how they see fit, e.g., vectorization */
36c7a4214aSPierre Jolivet       for (PetscInt k = 0; k < N; k++) ptr[j+k*M] = this->get_coef(rows[j],cols[k]);
37c7a4214aSPierre Jolivet     PetscFunctionReturnVoid();
38c7a4214aSPierre Jolivet   }
39c7a4214aSPierre Jolivet };
40c7a4214aSPierre Jolivet 
41c7a4214aSPierre Jolivet int main(int argc,char **argv)
42c7a4214aSPierre Jolivet {
43c7a4214aSPierre Jolivet   Mat            A,B,P,R;
44c7a4214aSPierre Jolivet   PetscInt       m = 100,dim = 3,M,begin = 0;
45c7a4214aSPierre Jolivet   PetscMPIInt    size;
46c7a4214aSPierre Jolivet   PetscReal      *coords,*gcoords,norm,epsilon,relative;
47c7a4214aSPierre Jolivet   PetscBool      sym = PETSC_FALSE;
48c7a4214aSPierre Jolivet   PetscRandom    rdm;
49c7a4214aSPierre Jolivet   MatHtoolKernel kernel = GenEntry;
50c7a4214aSPierre Jolivet   MyIMatrix      *imatrix;
51c7a4214aSPierre Jolivet   PetscErrorCode ierr;
52c7a4214aSPierre Jolivet 
53c7a4214aSPierre Jolivet   ierr = PetscInitialize(&argc,&argv,(char*)NULL,help);if (ierr) return ierr;
54c7a4214aSPierre Jolivet   ierr = PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL);CHKERRQ(ierr);
55c7a4214aSPierre Jolivet   ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr);
56c7a4214aSPierre Jolivet   ierr = PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL);CHKERRQ(ierr);
57c7a4214aSPierre Jolivet   ierr = PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL);CHKERRQ(ierr);
58c7a4214aSPierre Jolivet   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
59c7a4214aSPierre Jolivet   M = size*m;
60c7a4214aSPierre Jolivet   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
61c7a4214aSPierre Jolivet   ierr = PetscMalloc1(m*dim,&coords);CHKERRQ(ierr);
62c7a4214aSPierre Jolivet   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr);
63c7a4214aSPierre Jolivet   ierr = PetscRandomGetValuesReal(rdm,m*dim,coords);CHKERRQ(ierr);
64c7a4214aSPierre Jolivet   ierr = PetscCalloc1(M*dim,&gcoords);CHKERRQ(ierr);
65c7a4214aSPierre Jolivet   ierr = MPI_Exscan(&m,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
66c7a4214aSPierre Jolivet   ierr = PetscArraycpy(gcoords+begin*dim,coords,m*dim);CHKERRQ(ierr);
67c7a4214aSPierre Jolivet   ierr = MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
68c7a4214aSPierre Jolivet   imatrix = new MyIMatrix(M,M,dim,gcoords);
69c7a4214aSPierre Jolivet   ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,NULL,imatrix,&A);CHKERRQ(ierr); /* block-wise assembly using htool::IMatrix<PetscScalar>::copy_submatrix() */
70c7a4214aSPierre Jolivet   ierr = MatSetOption(A,MAT_SYMMETRIC,sym);CHKERRQ(ierr);
71c7a4214aSPierre Jolivet   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
72c7a4214aSPierre Jolivet   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
73c7a4214aSPierre Jolivet   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
74c7a4214aSPierre Jolivet   ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr);
75c7a4214aSPierre Jolivet   ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B);CHKERRQ(ierr); /* entry-wise assembly using GenEntry() */
76c7a4214aSPierre Jolivet   ierr = MatSetOption(B,MAT_SYMMETRIC,sym);CHKERRQ(ierr);
77c7a4214aSPierre Jolivet   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
78c7a4214aSPierre Jolivet   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79c7a4214aSPierre Jolivet   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80c7a4214aSPierre Jolivet   ierr = MatViewFromOptions(B,NULL,"-B_view");CHKERRQ(ierr);
81c7a4214aSPierre Jolivet   ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P);CHKERRQ(ierr);
82c7a4214aSPierre Jolivet   ierr = MatNorm(P,NORM_FROBENIUS,&relative);CHKERRQ(ierr);
83c7a4214aSPierre Jolivet   ierr = MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R);CHKERRQ(ierr);
84c7a4214aSPierre Jolivet   ierr = MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
85c7a4214aSPierre Jolivet   ierr = MatNorm(R,NORM_INFINITY,&norm);CHKERRQ(ierr);
86c7a4214aSPierre Jolivet   if (PetscAbsReal(norm/relative) > epsilon) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon);
87c7a4214aSPierre Jolivet   ierr = MatDestroy(&B);CHKERRQ(ierr);
88c7a4214aSPierre Jolivet   ierr = MatDestroy(&R);CHKERRQ(ierr);
89c7a4214aSPierre Jolivet   ierr = MatDestroy(&P);CHKERRQ(ierr);
90c7a4214aSPierre Jolivet   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
91c7a4214aSPierre Jolivet   ierr = MatDestroy(&A);CHKERRQ(ierr);
92c7a4214aSPierre Jolivet   ierr = PetscFree(gcoords);CHKERRQ(ierr);
93c7a4214aSPierre Jolivet   ierr = PetscFree(coords);CHKERRQ(ierr);
94c7a4214aSPierre Jolivet   delete imatrix;
95c7a4214aSPierre Jolivet   ierr = PetscFinalize();
96c7a4214aSPierre Jolivet   return ierr;
97c7a4214aSPierre Jolivet }
98c7a4214aSPierre Jolivet 
99c7a4214aSPierre Jolivet /*TEST
100c7a4214aSPierre Jolivet 
101c7a4214aSPierre Jolivet    build:
102c7a4214aSPierre Jolivet       requires: htool
103c7a4214aSPierre Jolivet    test:
104c7a4214aSPierre Jolivet       requires: htool
105c7a4214aSPierre Jolivet       suffix: 2
106c7a4214aSPierre Jolivet       nsize: 4
107c7a4214aSPierre Jolivet       args: -m_local 120 -mat_htool_epsilon 1.0e-2 -symmetric {{false true}shared output}
108c7a4214aSPierre Jolivet       output_file: output/ex101.out
109c7a4214aSPierre Jolivet 
110c7a4214aSPierre Jolivet TEST*/
111