Actual source code: power.c

  1: #include <petsc/private/matimpl.h>

  3: static PetscErrorCode MatColoringApply_Power(MatColoring mc, ISColoring *iscoloring)
  4: {
  5:   Mat         m = mc->mat, mp, ms;
  6:   MatColoring imc;
  7:   PetscInt    i;
  8:   const char *optionsprefix;

 10:   PetscFunctionBegin;
 11:   /* square the matrix repeatedly if necessary */
 12:   if (mc->dist == 1) {
 13:     mp = m;
 14:   } else {
 15:     PetscCall(MatMatMult(m, m, MAT_INITIAL_MATRIX, 2.0, &mp));
 16:     for (i = 2; i < mc->dist; i++) {
 17:       ms = mp;
 18:       PetscCall(MatMatMult(m, ms, MAT_INITIAL_MATRIX, 2.0, &mp));
 19:       PetscCall(MatDestroy(&ms));
 20:     }
 21:   }
 22:   PetscCall(MatColoringCreate(mp, &imc));
 23:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)mc, &optionsprefix));
 24:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)imc, optionsprefix));
 25:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)imc, "power_"));
 26:   PetscCall(MatColoringSetType(imc, MATCOLORINGGREEDY));
 27:   PetscCall(MatColoringSetDistance(imc, 1));
 28:   PetscCall(MatColoringSetWeightType(imc, mc->weight_type));
 29:   PetscCall(MatColoringSetFromOptions(imc));
 30:   PetscCall(MatColoringApply(imc, iscoloring));
 31:   PetscCall(MatColoringDestroy(&imc));
 32:   if (mp != m) PetscCall(MatDestroy(&mp));
 33:   PetscFunctionReturn(PETSC_SUCCESS);
 34: }

 36: /*MC
 37:   MATCOLORINGPOWER - Take the matrix's nth power, then do one-coloring on it.

 39:    Level: beginner

 41:    Notes:
 42:    This is merely a trivial test algorithm.

 44:    Supports any distance coloring.

 46: .seealso: `MatColoring`, `MatColoringType`, `MatColoringCreate()`, `MatColoring`, `MatColoringSetType()`
 47: M*/
 48: PETSC_EXTERN PetscErrorCode MatColoringCreate_Power(MatColoring mc)
 49: {
 50:   PetscFunctionBegin;
 51:   mc->ops->apply          = MatColoringApply_Power;
 52:   mc->ops->view           = NULL;
 53:   mc->ops->destroy        = NULL;
 54:   mc->ops->setfromoptions = NULL;
 55:   PetscFunctionReturn(PETSC_SUCCESS);
 56: }