Actual source code: packimpl.h

  1: #pragma once

  3: #include <petscdmcomposite.h>
  4: #include <petsc/private/dmimpl.h>

  6: /*
  7:    rstart is where an array/subvector starts in the global parallel vector, so arrays
  8:    rstarts are meaningless (and set to the previous one) except on the processor where the array lives
  9: */

 11: struct DMCompositeLink {
 12:   struct DMCompositeLink *next;
 13:   PetscInt                n;       /* number of owned */
 14:   PetscInt                rstart;  /* rstart is relative to this process */
 15:   PetscInt                grstart; /* grstart is relative to all processes */
 16:   PetscInt                nlocal;

 18:   /* only used for DMCOMPOSITE_DM */
 19:   PetscInt *grstarts; /* global row for first unknown of this DM on each process */
 20:   DM        dm;
 21: };

 23: typedef struct {
 24:   PetscInt                n, N, rstart; /* rstart is relative to all processors, n unknowns owned by this process, N is total unknowns */
 25:   PetscInt                nghost;       /* number of all local entries (includes DMDA ghost points) */
 26:   PetscInt                nDM, nmine;   /* how many DM's and separate redundant arrays used to build DM(nmine is ones on this process) */
 27:   PetscBool               setup;        /* after this is set, cannot add new links to the DM*/
 28:   struct DMCompositeLink *next;

 30:   PetscErrorCode (*FormCoupleLocations)(DM, Mat, PetscInt *, PetscInt *, PetscInt, PetscInt, PetscInt, PetscInt);
 31: } DM_Composite;

 33: PETSC_INTERN PetscErrorCode DMCreateMatrix_Composite(DM, Mat *);