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