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 549566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)NULL,help)); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL)); 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL)); 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL)); 599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 60c7a4214aSPierre Jolivet M = size*m; 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m*dim,&coords)); 639566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); 649566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValuesReal(rdm,m*dim,coords)); 659566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(M*dim,&gcoords)); 669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Exscan(&m,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 679566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(gcoords+begin*dim,coords,m*dim)); 681c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD)); 69c7a4214aSPierre Jolivet imatrix = new MyIMatrix(M,M,dim,gcoords); 709566063dSJacob Faibussowitsch PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,NULL,imatrix,&A)); /* block-wise assembly using htool::IMatrix<PetscScalar>::copy_submatrix() */ 719566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_SYMMETRIC,sym)); 729566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 739566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 749566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 759566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-A_view")); 769566063dSJacob Faibussowitsch PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B)); /* entry-wise assembly using GenEntries() */ 779566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_SYMMETRIC,sym)); 789566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 799566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 809566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 819566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(B,NULL,"-B_view")); 829566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P)); 839566063dSJacob Faibussowitsch PetscCall(MatNorm(P,NORM_FROBENIUS,&relative)); 849566063dSJacob Faibussowitsch PetscCall(MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R)); 859566063dSJacob Faibussowitsch PetscCall(MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN)); 869566063dSJacob Faibussowitsch PetscCall(MatNorm(R,NORM_INFINITY,&norm)); 87*08401ef6SPierre Jolivet PetscCheck(PetscAbsReal(norm/relative) <= epsilon,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon); 889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 919566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 939566063dSJacob Faibussowitsch PetscCall(PetscFree(gcoords)); 949566063dSJacob Faibussowitsch PetscCall(PetscFree(coords)); 95c7a4214aSPierre Jolivet delete imatrix; 969566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 97b122ec5aSJacob 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