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 } 21*cab00e0dSPierre Jolivet class MyIMatrix : public htool::VirtualGenerator<PetscScalar> { 22c7a4214aSPierre Jolivet private: 23c7a4214aSPierre Jolivet PetscReal *coords; 24c7a4214aSPierre Jolivet PetscInt sdim; 25c7a4214aSPierre Jolivet public: 26*cab00e0dSPierre 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 PetscErrorCode ierr; 54c7a4214aSPierre Jolivet 55c7a4214aSPierre Jolivet ierr = PetscInitialize(&argc,&argv,(char*)NULL,help);if (ierr) return ierr; 56c7a4214aSPierre Jolivet ierr = PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL);CHKERRQ(ierr); 57c7a4214aSPierre Jolivet ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr); 58c7a4214aSPierre Jolivet ierr = PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL);CHKERRQ(ierr); 59c7a4214aSPierre Jolivet ierr = PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL);CHKERRQ(ierr); 60c7a4214aSPierre Jolivet ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 61c7a4214aSPierre Jolivet M = size*m; 62c7a4214aSPierre Jolivet ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 63c7a4214aSPierre Jolivet ierr = PetscMalloc1(m*dim,&coords);CHKERRQ(ierr); 64c7a4214aSPierre Jolivet ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr); 65c7a4214aSPierre Jolivet ierr = PetscRandomGetValuesReal(rdm,m*dim,coords);CHKERRQ(ierr); 66c7a4214aSPierre Jolivet ierr = PetscCalloc1(M*dim,&gcoords);CHKERRQ(ierr); 67c7a4214aSPierre Jolivet ierr = MPI_Exscan(&m,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 68c7a4214aSPierre Jolivet ierr = PetscArraycpy(gcoords+begin*dim,coords,m*dim);CHKERRQ(ierr); 69c7a4214aSPierre Jolivet ierr = MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 70c7a4214aSPierre Jolivet imatrix = new MyIMatrix(M,M,dim,gcoords); 71c7a4214aSPierre 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() */ 72c7a4214aSPierre Jolivet ierr = MatSetOption(A,MAT_SYMMETRIC,sym);CHKERRQ(ierr); 73c7a4214aSPierre Jolivet ierr = MatSetFromOptions(A);CHKERRQ(ierr); 74c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 75c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 76c7a4214aSPierre Jolivet ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr); 7798e73e17SPierre Jolivet ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B);CHKERRQ(ierr); /* entry-wise assembly using GenEntries() */ 78c7a4214aSPierre Jolivet ierr = MatSetOption(B,MAT_SYMMETRIC,sym);CHKERRQ(ierr); 79c7a4214aSPierre Jolivet ierr = MatSetFromOptions(B);CHKERRQ(ierr); 80c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 81c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 82c7a4214aSPierre Jolivet ierr = MatViewFromOptions(B,NULL,"-B_view");CHKERRQ(ierr); 83c7a4214aSPierre Jolivet ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P);CHKERRQ(ierr); 84c7a4214aSPierre Jolivet ierr = MatNorm(P,NORM_FROBENIUS,&relative);CHKERRQ(ierr); 85c7a4214aSPierre Jolivet ierr = MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R);CHKERRQ(ierr); 86c7a4214aSPierre Jolivet ierr = MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 87c7a4214aSPierre Jolivet ierr = MatNorm(R,NORM_INFINITY,&norm);CHKERRQ(ierr); 88c7a4214aSPierre 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); 89c7a4214aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 90c7a4214aSPierre Jolivet ierr = MatDestroy(&R);CHKERRQ(ierr); 91c7a4214aSPierre Jolivet ierr = MatDestroy(&P);CHKERRQ(ierr); 92c7a4214aSPierre Jolivet ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 93c7a4214aSPierre Jolivet ierr = MatDestroy(&A);CHKERRQ(ierr); 94c7a4214aSPierre Jolivet ierr = PetscFree(gcoords);CHKERRQ(ierr); 95c7a4214aSPierre Jolivet ierr = PetscFree(coords);CHKERRQ(ierr); 96c7a4214aSPierre Jolivet delete imatrix; 97c7a4214aSPierre Jolivet ierr = PetscFinalize(); 98c7a4214aSPierre Jolivet return ierr; 99c7a4214aSPierre Jolivet } 100c7a4214aSPierre Jolivet 101c7a4214aSPierre Jolivet /*TEST 102c7a4214aSPierre Jolivet 103c7a4214aSPierre Jolivet build: 104c7a4214aSPierre Jolivet requires: htool 105c7a4214aSPierre Jolivet test: 106c7a4214aSPierre Jolivet requires: htool 107c7a4214aSPierre Jolivet suffix: 2 108c7a4214aSPierre Jolivet nsize: 4 109c7a4214aSPierre Jolivet args: -m_local 120 -mat_htool_epsilon 1.0e-2 -symmetric {{false true}shared output} 110c7a4214aSPierre Jolivet output_file: output/ex101.out 111c7a4214aSPierre Jolivet 112c7a4214aSPierre Jolivet TEST*/ 113