xref: /petsc/src/mat/tests/ex268.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
11dd4f53aSPierre Jolivet static char help[] = "Tests MATFACTORHTOOL\n\n";
21dd4f53aSPierre Jolivet 
31dd4f53aSPierre Jolivet #include <petscmat.h>
41dd4f53aSPierre Jolivet 
GenEntries(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt * J,const PetscInt * K,PetscScalar * ptr,PetscCtx ctx)5*2a8381b2SBarry Smith static PetscErrorCode GenEntries(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *J, const PetscInt *K, PetscScalar *ptr, PetscCtx ctx)
61dd4f53aSPierre Jolivet {
71dd4f53aSPierre Jolivet   PetscInt  d, j, k;
81dd4f53aSPierre Jolivet   PetscReal diff = 0.0, *coords = (PetscReal *)(ctx);
91dd4f53aSPierre Jolivet 
101dd4f53aSPierre Jolivet   PetscFunctionBeginUser;
111dd4f53aSPierre Jolivet   for (j = 0; j < M; j++) {
121dd4f53aSPierre Jolivet     for (k = 0; k < N; k++) {
131dd4f53aSPierre Jolivet       diff = 0.0;
141dd4f53aSPierre 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]);
151dd4f53aSPierre Jolivet       ptr[j + M * k] = 1.0 / (1.0e-1 + PetscSqrtReal(diff));
161dd4f53aSPierre Jolivet     }
171dd4f53aSPierre Jolivet   }
181dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
191dd4f53aSPierre Jolivet }
201dd4f53aSPierre Jolivet 
main(int argc,char ** argv)211dd4f53aSPierre Jolivet int main(int argc, char **argv)
221dd4f53aSPierre Jolivet {
231dd4f53aSPierre Jolivet   Mat               A, Ad, F, Fd, X, Xd, B;
241dd4f53aSPierre Jolivet   Vec               x, xd, b;
251dd4f53aSPierre Jolivet   PetscInt          m = 100, dim = 3, M, K = 10, begin, n = 0;
261dd4f53aSPierre Jolivet   PetscMPIInt       size;
271dd4f53aSPierre Jolivet   PetscReal        *coords, *gcoords, norm, epsilon;
281dd4f53aSPierre Jolivet   MatHtoolKernelFn *kernel = GenEntries;
291dd4f53aSPierre Jolivet   PetscBool         flg, sym = PETSC_FALSE;
301dd4f53aSPierre Jolivet   PetscRandom       rdm;
311dd4f53aSPierre Jolivet   MatSolverType     type;
321dd4f53aSPierre Jolivet 
331dd4f53aSPierre Jolivet   PetscFunctionBeginUser;
341dd4f53aSPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, (char *)NULL, help));
351dd4f53aSPierre Jolivet   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m_local", &m, NULL));
361dd4f53aSPierre Jolivet   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n_local", &n, NULL));
371dd4f53aSPierre Jolivet   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
381dd4f53aSPierre Jolivet   PetscCall(PetscOptionsGetInt(NULL, NULL, "-K", &K, NULL));
391dd4f53aSPierre Jolivet   PetscCall(PetscOptionsGetBool(NULL, NULL, "-symmetric", &sym, NULL));
401dd4f53aSPierre Jolivet   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mat_htool_epsilon", &epsilon, NULL));
411dd4f53aSPierre Jolivet   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
421dd4f53aSPierre Jolivet   M = size * m;
431dd4f53aSPierre Jolivet   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
441dd4f53aSPierre Jolivet   PetscCall(PetscMalloc1(m * dim, &coords));
451dd4f53aSPierre Jolivet   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
461dd4f53aSPierre Jolivet   PetscCall(PetscRandomGetValuesReal(rdm, m * dim, coords));
471dd4f53aSPierre Jolivet   PetscCall(PetscCalloc1(M * dim, &gcoords));
481dd4f53aSPierre Jolivet   PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, PETSC_DECIDE, M, K, NULL, &B));
491dd4f53aSPierre Jolivet   PetscCall(MatSetRandom(B, rdm));
501dd4f53aSPierre Jolivet   PetscCall(MatGetOwnershipRange(B, &begin, NULL));
511dd4f53aSPierre Jolivet   PetscCall(PetscArraycpy(gcoords + begin * dim, coords, m * dim));
521dd4f53aSPierre Jolivet   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, gcoords, M * dim, MPIU_REAL, MPI_SUM, PETSC_COMM_WORLD));
531dd4f53aSPierre Jolivet   PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD, m, m, M, M, dim, coords, coords, kernel, gcoords, &A));
541dd4f53aSPierre Jolivet   PetscCall(MatSetOption(A, MAT_SYMMETRIC, sym));
551dd4f53aSPierre Jolivet   PetscCall(MatSetFromOptions(A));
561dd4f53aSPierre Jolivet   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
571dd4f53aSPierre Jolivet   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
581dd4f53aSPierre Jolivet   PetscCall(MatConvert(A, MATDENSE, MAT_INITIAL_MATRIX, &Ad));
591dd4f53aSPierre Jolivet   PetscCall(MatMultEqual(A, Ad, 10, &flg));
601dd4f53aSPierre Jolivet   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Ax != Adx");
611dd4f53aSPierre Jolivet   PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, PETSC_DECIDE, M, K, NULL, &X));
621dd4f53aSPierre Jolivet   PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, PETSC_DECIDE, M, K, NULL, &Xd));
631dd4f53aSPierre Jolivet   PetscCall(MatViewFromOptions(A, NULL, "-A"));
641dd4f53aSPierre Jolivet   PetscCall(MatViewFromOptions(Ad, NULL, "-Ad"));
651dd4f53aSPierre Jolivet   PetscCall(MatViewFromOptions(B, NULL, "-B"));
661dd4f53aSPierre Jolivet   for (PetscInt i = 0; i < 2; ++i) {
671dd4f53aSPierre Jolivet     PetscCall(MatGetFactor(A, MATSOLVERHTOOL, i == 0 ? MAT_FACTOR_LU : MAT_FACTOR_CHOLESKY, &F));
681dd4f53aSPierre Jolivet     PetscCall(MatGetFactor(Ad, MATSOLVERPETSC, i == 0 ? MAT_FACTOR_LU : MAT_FACTOR_CHOLESKY, &Fd));
691dd4f53aSPierre Jolivet     PetscCall(MatFactorGetSolverType(F, &type));
701dd4f53aSPierre Jolivet     PetscCall(PetscStrncmp(type, MATSOLVERHTOOL, 5, &flg));
711dd4f53aSPierre Jolivet     PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "MATSOLVERHTOOL != htool");
721dd4f53aSPierre Jolivet     if (i == 0) {
731dd4f53aSPierre Jolivet       PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL));
741dd4f53aSPierre Jolivet       PetscCall(MatLUFactorNumeric(F, A, NULL));
751dd4f53aSPierre Jolivet       PetscCall(MatLUFactorSymbolic(Fd, Ad, NULL, NULL, NULL));
761dd4f53aSPierre Jolivet       PetscCall(MatLUFactorNumeric(Fd, Ad, NULL));
771dd4f53aSPierre Jolivet     } else {
781dd4f53aSPierre Jolivet       PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL));
791dd4f53aSPierre Jolivet       PetscCall(MatCholeskyFactorNumeric(F, A, NULL));
801dd4f53aSPierre Jolivet       PetscCall(MatCholeskyFactorSymbolic(Fd, Ad, NULL, NULL));
811dd4f53aSPierre Jolivet       PetscCall(MatCholeskyFactorNumeric(Fd, Ad, NULL));
821dd4f53aSPierre Jolivet     }
831dd4f53aSPierre Jolivet     PetscCall(MatMatSolve(F, B, X));
841dd4f53aSPierre Jolivet     PetscCall(MatMatSolve(Fd, B, Xd));
851dd4f53aSPierre Jolivet     PetscCall(MatViewFromOptions(X, NULL, "-X"));
861dd4f53aSPierre Jolivet     PetscCall(MatViewFromOptions(Xd, NULL, "-Xd"));
871dd4f53aSPierre Jolivet     PetscCall(MatAXPY(Xd, -1.0, X, SAME_NONZERO_PATTERN));
881dd4f53aSPierre Jolivet     PetscCall(MatNorm(Xd, NORM_INFINITY, &norm));
891dd4f53aSPierre Jolivet     PetscCall(MatViewFromOptions(Xd, NULL, "-MatMatSolve"));
901dd4f53aSPierre Jolivet     if (norm > 0.01) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: norm of residual for MatMatSolve %g\n", (double)norm));
911dd4f53aSPierre Jolivet     if (!PetscDefined(USE_COMPLEX) || i == 0) {
921dd4f53aSPierre Jolivet       PetscCall(MatMatSolveTranspose(F, B, X));
931dd4f53aSPierre Jolivet       PetscCall(MatMatSolveTranspose(Fd, B, Xd));
941dd4f53aSPierre Jolivet       PetscCall(MatAXPY(Xd, -1.0, X, SAME_NONZERO_PATTERN));
951dd4f53aSPierre Jolivet       PetscCall(MatNorm(Xd, NORM_INFINITY, &norm));
961dd4f53aSPierre Jolivet       PetscCall(MatViewFromOptions(Xd, NULL, "-MatMatSolveTranspose"));
971dd4f53aSPierre Jolivet       if (norm > 0.01) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: norm of residual for MatMatSolveTranspose %g\n", (double)norm));
981dd4f53aSPierre Jolivet     }
991dd4f53aSPierre Jolivet     PetscCall(MatDenseGetColumnVecRead(B, 0, &b));
1001dd4f53aSPierre Jolivet     PetscCall(MatDenseGetColumnVecWrite(X, 0, &x));
1011dd4f53aSPierre Jolivet     PetscCall(MatDenseGetColumnVecWrite(Xd, 0, &xd));
1021dd4f53aSPierre Jolivet     PetscCall(MatSolve(F, b, x));
1031dd4f53aSPierre Jolivet     PetscCall(MatSolve(Fd, b, xd));
1041dd4f53aSPierre Jolivet     PetscCall(VecAXPY(xd, -1.0, x));
1051dd4f53aSPierre Jolivet     PetscCall(VecNorm(xd, NORM_INFINITY, &norm));
1061dd4f53aSPierre Jolivet     PetscCall(MatViewFromOptions(Xd, NULL, "-MatSolve"));
1071dd4f53aSPierre Jolivet     if (norm > 0.01) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: norm of residual for MatSolve %g\n", (double)norm));
1081dd4f53aSPierre Jolivet     if (!PetscDefined(USE_COMPLEX) || i == 0) {
1091dd4f53aSPierre Jolivet       PetscCall(MatSolveTranspose(F, b, x));
1101dd4f53aSPierre Jolivet       PetscCall(MatSolveTranspose(Fd, b, xd));
1111dd4f53aSPierre Jolivet       PetscCall(VecAXPY(xd, -1.0, x));
1121dd4f53aSPierre Jolivet       PetscCall(VecNorm(xd, NORM_INFINITY, &norm));
1131dd4f53aSPierre Jolivet       PetscCall(MatViewFromOptions(Xd, NULL, "-MatSolveTranspose"));
1141dd4f53aSPierre Jolivet       if (norm > 0.01) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: norm of residual for MatSolveTranspose %g\n", (double)norm));
1151dd4f53aSPierre Jolivet     }
1161dd4f53aSPierre Jolivet     PetscCall(MatDenseRestoreColumnVecWrite(Xd, 0, &xd));
1171dd4f53aSPierre Jolivet     PetscCall(MatDenseRestoreColumnVecWrite(X, 0, &x));
1181dd4f53aSPierre Jolivet     PetscCall(MatDenseRestoreColumnVecRead(B, 0, &b));
1191dd4f53aSPierre Jolivet     PetscCall(MatDestroy(&Fd));
1201dd4f53aSPierre Jolivet     PetscCall(MatDestroy(&F));
1211dd4f53aSPierre Jolivet   }
1221dd4f53aSPierre Jolivet   PetscCall(MatDestroy(&Xd));
1231dd4f53aSPierre Jolivet   PetscCall(MatDestroy(&X));
1241dd4f53aSPierre Jolivet   PetscCall(PetscRandomDestroy(&rdm));
1251dd4f53aSPierre Jolivet   PetscCall(MatDestroy(&Ad));
1261dd4f53aSPierre Jolivet   PetscCall(MatDestroy(&A));
1271dd4f53aSPierre Jolivet   PetscCall(MatDestroy(&B));
1281dd4f53aSPierre Jolivet   PetscCall(PetscFree(gcoords));
1291dd4f53aSPierre Jolivet   PetscCall(PetscFree(coords));
1301dd4f53aSPierre Jolivet   PetscCall(PetscFinalize());
1311dd4f53aSPierre Jolivet   return 0;
1321dd4f53aSPierre Jolivet }
1331dd4f53aSPierre Jolivet 
1341dd4f53aSPierre Jolivet /*TEST
1351dd4f53aSPierre Jolivet 
1361dd4f53aSPierre Jolivet    build:
1371dd4f53aSPierre Jolivet       requires: htool
1381dd4f53aSPierre Jolivet 
1391dd4f53aSPierre Jolivet    test:
1401dd4f53aSPierre Jolivet       requires: htool
1411dd4f53aSPierre Jolivet       suffix: 1
1421dd4f53aSPierre Jolivet       nsize: 1
1431dd4f53aSPierre Jolivet       args: -mat_htool_epsilon 1.0e-11
1441dd4f53aSPierre Jolivet       output_file: output/empty.out
1451dd4f53aSPierre Jolivet 
1461dd4f53aSPierre Jolivet TEST*/
147