xref: /petsc/src/mat/tests/ex246.cxx (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
698e73e17SPierre Jolivet static PetscErrorCode GenEntries(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *J,const PetscInt *K,PetscScalar *ptr,void *ctx)
7c7a4214aSPierre Jolivet {
898e73e17SPierre Jolivet   PetscInt  d,j,k;
9c7a4214aSPierre Jolivet   PetscReal diff = 0.0,*coords = (PetscReal*)(ctx);
10c7a4214aSPierre Jolivet 
11c7a4214aSPierre Jolivet   PetscFunctionBeginUser;
1298e73e17SPierre Jolivet   for (j = 0; j < M; j++) {
1398e73e17SPierre Jolivet     for (k = 0; k < N; k++) {
1498e73e17SPierre Jolivet       diff = 0.0;
1598e73e17SPierre Jolivet       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]);
1698e73e17SPierre Jolivet       ptr[j+M*k] = 1.0/(1.0e-2 + PetscSqrtReal(diff));
17c7a4214aSPierre Jolivet     }
1898e73e17SPierre Jolivet   }
1998e73e17SPierre Jolivet   PetscFunctionReturn(0);
2098e73e17SPierre Jolivet }
21cab00e0dSPierre Jolivet class MyIMatrix : public htool::VirtualGenerator<PetscScalar> {
22c7a4214aSPierre Jolivet   private:
23c7a4214aSPierre Jolivet   PetscReal *coords;
24c7a4214aSPierre Jolivet   PetscInt  sdim;
25c7a4214aSPierre Jolivet   public:
26cab00e0dSPierre Jolivet   MyIMatrix(PetscInt M,PetscInt N,PetscInt spacedim,PetscReal* gcoords) : htool::VirtualGenerator<PetscScalar>(M,N),coords(gcoords),sdim(spacedim) { }
27c7a4214aSPierre Jolivet 
2898e73e17SPierre Jolivet   void copy_submatrix(PetscInt M,PetscInt N,const PetscInt *J,const PetscInt *K,PetscScalar *ptr) const override
2970a7d78aSStefano Zampini   {
30c7a4214aSPierre Jolivet     PetscReal diff = 0.0;
31c7a4214aSPierre Jolivet 
32c7a4214aSPierre Jolivet     PetscFunctionBeginUser;
33c7a4214aSPierre Jolivet     for (PetscInt j = 0; j < M; j++) /* could be optimized by the user how they see fit, e.g., vectorization */
3498e73e17SPierre Jolivet       for (PetscInt k = 0; k < N; k++) {
3598e73e17SPierre Jolivet         diff = 0.0;
3698e73e17SPierre Jolivet         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]);
3798e73e17SPierre Jolivet         ptr[j+M*k] = 1.0/(1.0e-2 + PetscSqrtReal(diff));
3898e73e17SPierre Jolivet       }
39c7a4214aSPierre Jolivet     PetscFunctionReturnVoid();
40c7a4214aSPierre Jolivet   }
41c7a4214aSPierre Jolivet };
42c7a4214aSPierre Jolivet 
43c7a4214aSPierre Jolivet int main(int argc,char **argv)
44c7a4214aSPierre Jolivet {
45c7a4214aSPierre Jolivet   Mat            A,B,P,R;
46c7a4214aSPierre Jolivet   PetscInt       m = 100,dim = 3,M,begin = 0;
47c7a4214aSPierre Jolivet   PetscMPIInt    size;
48c7a4214aSPierre Jolivet   PetscReal      *coords,*gcoords,norm,epsilon,relative;
49c7a4214aSPierre Jolivet   PetscBool      sym = PETSC_FALSE;
50c7a4214aSPierre Jolivet   PetscRandom    rdm;
5198e73e17SPierre Jolivet   MatHtoolKernel kernel = GenEntries;
52c7a4214aSPierre Jolivet   MyIMatrix      *imatrix;
53c7a4214aSPierre Jolivet 
54*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)NULL,help));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL));
595f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
60c7a4214aSPierre Jolivet   M = size*m;
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(m*dim,&coords));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rdm));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomGetValuesReal(rdm,m*dim,coords));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(M*dim,&gcoords));
665f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Exscan(&m,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(gcoords+begin*dim,coords,m*dim));
685f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD));
69c7a4214aSPierre Jolivet   imatrix = new MyIMatrix(M,M,dim,gcoords);
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,NULL,imatrix,&A)); /* block-wise assembly using htool::IMatrix<PetscScalar>::copy_submatrix() */
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,sym));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(A,NULL,"-A_view"));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B)); /* entry-wise assembly using GenEntries() */
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(B,MAT_SYMMETRIC,sym));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(B));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(B,NULL,"-B_view"));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(P,NORM_FROBENIUS,&relative));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(R,NORM_INFINITY,&norm));
872c71b3e2SJacob Faibussowitsch   PetscCheckFalse(PetscAbsReal(norm/relative) > epsilon,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon);
885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&R));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&P));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(gcoords));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(coords));
95c7a4214aSPierre Jolivet   delete imatrix;
96*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
97*b122ec5aSJacob Faibussowitsch   return 0;
98c7a4214aSPierre Jolivet }
99c7a4214aSPierre Jolivet 
100c7a4214aSPierre Jolivet /*TEST
101c7a4214aSPierre Jolivet 
102c7a4214aSPierre Jolivet    build:
103c7a4214aSPierre Jolivet       requires: htool
104c7a4214aSPierre Jolivet    test:
105c7a4214aSPierre Jolivet       requires: htool
106c7a4214aSPierre Jolivet       suffix: 2
107c7a4214aSPierre Jolivet       nsize: 4
108c7a4214aSPierre Jolivet       args: -m_local 120 -mat_htool_epsilon 1.0e-2 -symmetric {{false true}shared output}
109c7a4214aSPierre Jolivet       output_file: output/ex101.out
110c7a4214aSPierre Jolivet 
111c7a4214aSPierre Jolivet TEST*/
112