Actual source code: pcpatchimpl.h

  1: /*
  2:       Data structure used for Patch preconditioner.
  3: */
  4: #pragma once
  5: #include <petsc/private/pcimpl.h>
  6: #include <petsc/private/hashseti.h>
  7: #include <petsc/private/hashmapi.h>
  8: #include <petscksp.h>

 10: PETSC_INTERN PetscBool  PCPatchcite;
 11: PETSC_INTERN const char PCPatchCitation[];

 13: typedef struct {
 14:   /* Topology */
 15:   PCPatchConstructType ctype;                                           /* Algorithm for patch construction */
 16:   PetscErrorCode (*patchconstructop)(void *, DM, PetscInt, PetscHSetI); /* patch construction */
 17:   PetscErrorCode (*userpatchconstructionop)(PC, PetscInt *, IS **, IS *, PetscCtx ctx);
 18:   void        *userpatchconstructctx;
 19:   IS          *userIS;
 20:   PetscInt     npatch;       /* Number of patches */
 21:   PetscBool    user_patches; /* Flag for user construction of patches */
 22:   PetscInt     dim, codim;   /* Dimension or codimension of mesh points to loop over; only one of them can be set */
 23:   PetscSection cellCounts;   /* Maps patch -> # cells in patch */
 24:   IS           cells;        /* [patch][cell in patch]: Cell number */
 25:   IS           extFacets;
 26:   IS           intFacets;
 27:   IS           intFacetsToPatchCell; /* Support of interior facet in local patch point numbering: AKA which two cells touch the facet (in patch local numbering of cells) */
 28:   IS           extFacetsToPatchCell; /* Support of exterior facet in local patch point numbering: the one cell that touches the facet (in patch local numbering of cells) */
 29:   PetscSection intFacetCounts;
 30:   PetscSection extFacetCounts;
 31:   PetscSection cellNumbering; /* Plex: NULL Firedrake: Numbering of cells in DM */
 32:   PetscSection pointCounts;   /* Maps patch -> # points with dofs in patch */
 33:   IS           points;        /* [patch][point in patch]: Point number */
 34:   /* Dof layout */
 35:   PetscBool     combined;        /* Use a combined space with all fields */
 36:   PetscInt      nsubspaces;      /* Number of fields */
 37:   PetscSF       sectionSF;       /* Combined SF mapping process local to global */
 38:   PetscSection *dofSection;      /* ?? For each field, patch -> # dofs in patch */
 39:   PetscInt     *subspaceOffsets; /* Plex: NULL Firedrake: offset of each field in concatenated process local numbering for mixed spaces */
 40:   PetscInt    **cellNodeMap;     /* [field][cell][dof in cell]: global dofs in cell TODO Free this after its use in PCPatchCreateCellPatchDiscretisationInfo() */
 41:   IS            dofs;            /* [patch][cell in patch][dof in cell]: patch local dof */
 42:   IS            offs;            /* [patch][point in patch]: patch local offset (same layout as 'points', used for filling up patchSection) */
 43:   IS            dofsWithArtificial;
 44:   IS            offsWithArtificial;
 45:   IS            dofsWithAll;
 46:   IS            offsWithAll;
 47:   PetscSection  patchSection;             /* Maps points -> patch local dofs */
 48:   IS            globalBcNodes;            /* Global dofs constrained by global Dirichlet conditions TODO Replace these with process local constrained dofs */
 49:   IS            ghostBcNodes;             /* Global dofs constrained by global Dirichlet conditions on this process and possibly others (patch overlaps boundary) */
 50:   PetscSection  gtolCounts;               /* ?? Indices to extract from local to patch vectors */
 51:   PetscSection  gtolCountsWithArtificial; /* ?? Indices to extract from local to patch vectors including those with artificial bcs*/
 52:   PetscSection  gtolCountsWithAll;        /* ?? Indices to extract from local to patch vectors including those in artificial or global bcs*/
 53:   IS            gtol;
 54:   IS            gtolWithArtificial;
 55:   IS            gtolWithAll;
 56:   PetscInt     *bs;                   /* [field] block size per field (can come from global operators?) */
 57:   PetscInt     *nodesPerCell;         /* [field] Dofs per cell TODO Change "node" to "dof" everywhere */
 58:   PetscInt      totalDofsPerCell;     /* Dofs per cell counting all fields */
 59:   PetscHSetI    subspaces_to_exclude; /* If you don't want any other dofs from a particular subspace you can exclude them with this.
 60:                                                 Used for Vanka in Stokes, for example, to eliminate all pressure dofs not on the vertex
 61:                                                 you're building the patch around */
 62:   PetscInt      vankadim;             /* In Vanka construction, should we eliminate any entities of a certain dimension on the initial patch? */
 63:   PetscInt      ignoredim;            /* In Vanka construction, should we eliminate any entities of a certain dimension on the boundary? */
 64:   PetscInt      pardecomp_overlap;    /* In parallel decomposition construction, how much overlap? */
 65:   /* Patch system assembly */
 66:   PetscErrorCode (*usercomputeop)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *);
 67:   void *usercomputeopctx;
 68:   PetscErrorCode (*usercomputef)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *);
 69:   void *usercomputefctx;
 70:   /* Interior facet integrals: Jacobian */
 71:   PetscErrorCode (*usercomputeopintfacet)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *);
 72:   void *usercomputeopintfacetctx;
 73:   /* Interior facet integrals: Residual */
 74:   PetscErrorCode (*usercomputefintfacet)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *);
 75:   void *usercomputefintfacetctx;
 76:   /* Exterior facet integrals: Jacobian */
 77:   PetscErrorCode (*usercomputeopextfacet)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *);
 78:   void *usercomputeopextfacetctx;
 79:   /* Exterior facet integrals: Residual */
 80:   PetscErrorCode (*usercomputefextfacet)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *);
 81:   void           *usercomputefextfacetctx;
 82:   IS              cellIS;                   /* Temporary IS for each cell patch */
 83:   PetscBool       save_operators;           /* Save all operators (or create/destroy one at a time?) */
 84:   PetscBool       precomputeElementTensors; /* Precompute all element tensors (each cell is assembled exactly once)? */
 85:   IS              allCells;                 /* Unique cells in union of all patches */
 86:   IS              allIntFacets;             /* Unique interior facets in union of all patches */
 87:   PetscBool       partition_of_unity;       /* Weight updates by dof multiplicity? */
 88:   PetscBool       multiplicative;           /* Gauss-Seidel instead of Jacobi?  */
 89:   PCCompositeType local_composition_type;   /* locally additive or multiplicative? */
 90:   /* Patch solves */
 91:   Vec          cellMats;                           /* Cell element tensors */
 92:   PetscInt    *precomputedTensorLocations;         /* Locations of the precomputed tensors for each cell. */
 93:   Vec          intFacetMats;                       /* interior facet element tensors */
 94:   PetscInt    *precomputedIntFacetTensorLocations; /* Locations of the precomputed tensors for each interior facet. */
 95:   Mat         *mat;                                /* System matrix for each patch */
 96:   Mat         *matWithArtificial;                  /* System matrix including dofs with artificial bcs for each patch */
 97:   MatType      sub_mat_type;                       /* Matrix type for patch systems */
 98:   Vec          patchRHS, patchUpdate;              /* Work vectors for RHS and solution on each patch */
 99:   IS          *dofMappingWithoutToWithArtificial;
100:   IS          *dofMappingWithoutToWithAll;
101:   Vec          patchRHSWithArtificial;         /* like patchRHS but extra entries to include dofs with artificial bcs*/
102:   Vec         *patch_dof_weights;              /* Weighting for dof in each patch */
103:   Vec          localRHS, localUpdate;          /* ??? */
104:   Vec          dof_weights;                    /* In how many patches does each dof lie? */
105:   PetscBool    symmetrise_sweep;               /* Should we sweep forwards->backwards, backwards->forwards? */
106:   PetscBool    optionsSet;                     /* SetFromOptions was called on this PC */
107:   IS           iterationSet;                   /* Index set specifying how we iterate over patches */
108:   PetscInt     currentPatch;                   /* The current patch number when iterating */
109:   PetscObject *solver;                         /* Solvers for each patch TODO Do we need a new KSP for each patch? */
110:   PetscBool    denseinverse;                   /* Should the patch inverse by applied by computing the inverse and a matmult? (Skips KSP/PC etc...) */
111:   PetscErrorCode (*densesolve)(Mat, Vec, Vec); /* Matmult for dense solve (used with denseinverse) */
112:   PetscErrorCode (*setupsolver)(PC);
113:   PetscErrorCode (*applysolver)(PC, PetscInt, Vec, Vec);
114:   PetscErrorCode (*resetsolver)(PC);
115:   PetscErrorCode (*destroysolver)(PC);
116:   PetscErrorCode (*updatemultiplicative)(PC, PetscInt, PetscInt);
117:   /* Monitoring */
118:   PetscBool         viewPatches;     /* View information about patch construction */
119:   PetscBool         viewCells;       /* View cells for each patch */
120:   PetscViewer       viewerCells;     /*   Viewer for patch cells */
121:   PetscViewerFormat formatCells;     /*   Format for patch cells */
122:   PetscBool         viewIntFacets;   /* View intFacets for each patch */
123:   PetscViewer       viewerIntFacets; /*   Viewer for patch intFacets */
124:   PetscViewerFormat formatIntFacets; /*   Format for patch intFacets */
125:   PetscBool         viewExtFacets;   /* View extFacets for each patch */
126:   PetscViewer       viewerExtFacets; /*   Viewer for patch extFacets */
127:   PetscViewerFormat formatExtFacets; /*   Format for patch extFacets */
128:   PetscBool         viewPoints;      /* View points for each patch */
129:   PetscViewer       viewerPoints;    /*   Viewer for patch points */
130:   PetscViewerFormat formatPoints;    /*   Format for patch points */
131:   PetscBool         viewSection;     /* View global section for each patch */
132:   PetscViewer       viewerSection;   /*   Viewer for patch sections */
133:   PetscViewerFormat formatSection;   /*   Format for patch sections */
134:   PetscBool         viewMatrix;      /* View matrix for each patch */
135:   PetscViewer       viewerMatrix;    /*   Viewer for patch matrix */
136:   PetscViewerFormat formatMatrix;    /*   Format for patch matrix */
137:   /* Extra variables for SNESPATCH */
138:   Vec         patchState;        /* State vectors for patch solvers */
139:   Vec         patchStateWithAll; /* State vectors for patch solvers with all boundary data */
140:   Vec         localState;        /* Scatter vector for state */
141:   Vec         patchResidual;     /* Work vectors for patch residual evaluation*/
142:   const char *classname;         /* "snes" or "pc" for options */
143:   PetscBool   isNonlinear;       /* we need to do some things differently in nonlinear mode */
144: } PC_PATCH;

146: PETSC_EXTERN PetscLogEvent PC_Patch_CreatePatches;
147: PETSC_EXTERN PetscLogEvent PC_Patch_ComputeOp;
148: PETSC_EXTERN PetscLogEvent PC_Patch_Solve;
149: PETSC_EXTERN PetscLogEvent PC_Patch_Apply;
150: PETSC_EXTERN PetscLogEvent PC_Patch_Prealloc;

152: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode PCPatchComputeFunction_Internal(PC, Vec, Vec, PetscInt);
153: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode PCPatchComputeOperator_Internal(PC, Vec, Mat, PetscInt, PetscBool);
154: typedef enum {
155:   SCATTER_INTERIOR,
156:   SCATTER_WITHARTIFICIAL,
157:   SCATTER_WITHALL
158: } PatchScatterType;
159: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode PCPatch_ScatterLocal_Private(PC, PetscInt, Vec, Vec, InsertMode, ScatterMode, PatchScatterType);