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