xref: /petsc/src/mat/tests/ex246.cxx (revision 7d5fd1e4d9337468ad3f05b65b7facdcd2dfd2a4)
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   ierr = PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL);CHKERRQ(ierr);
57   ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr);
58   ierr = PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL);CHKERRQ(ierr);
59   ierr = PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL);CHKERRQ(ierr);
60   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
61   M = size*m;
62   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
63   ierr = PetscMalloc1(m*dim,&coords);CHKERRQ(ierr);
64   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr);
65   ierr = PetscRandomGetValuesReal(rdm,m*dim,coords);CHKERRQ(ierr);
66   ierr = PetscCalloc1(M*dim,&gcoords);CHKERRQ(ierr);
67   ierr = MPI_Exscan(&m,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
68   ierr = PetscArraycpy(gcoords+begin*dim,coords,m*dim);CHKERRQ(ierr);
69   ierr = MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
70   imatrix = new MyIMatrix(M,M,dim,gcoords);
71   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() */
72   ierr = MatSetOption(A,MAT_SYMMETRIC,sym);CHKERRQ(ierr);
73   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
74   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
75   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76   ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr);
77   ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B);CHKERRQ(ierr); /* entry-wise assembly using GenEntries() */
78   ierr = MatSetOption(B,MAT_SYMMETRIC,sym);CHKERRQ(ierr);
79   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
80   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
81   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
82   ierr = MatViewFromOptions(B,NULL,"-B_view");CHKERRQ(ierr);
83   ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P);CHKERRQ(ierr);
84   ierr = MatNorm(P,NORM_FROBENIUS,&relative);CHKERRQ(ierr);
85   ierr = MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R);CHKERRQ(ierr);
86   ierr = MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
87   ierr = MatNorm(R,NORM_INFINITY,&norm);CHKERRQ(ierr);
88   if (PetscAbsReal(norm/relative) > epsilon) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon);
89   ierr = MatDestroy(&B);CHKERRQ(ierr);
90   ierr = MatDestroy(&R);CHKERRQ(ierr);
91   ierr = MatDestroy(&P);CHKERRQ(ierr);
92   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
93   ierr = MatDestroy(&A);CHKERRQ(ierr);
94   ierr = PetscFree(gcoords);CHKERRQ(ierr);
95   ierr = PetscFree(coords);CHKERRQ(ierr);
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