Actual source code: pcis.h

petsc-3.12.5 2020-03-29
Report Typos and Errors


  5:  #include <petsc/private/pcimpl.h>
  6:  #include <../src/mat/impls/is/matis.h>
  7:  #include <petscksp.h>

  9: /*
 10:    Context (data structure) common for all Iterative Substructuring preconditioners.
 11: */

 13: typedef struct {

 15:   /* In naming the variables, we adopted the following convention: */
 16:   /* * B - stands for interface nodes;                             */
 17:   /* * I - stands for interior nodes;                              */
 18:   /* * D - stands for Dirichlet (by extension, refers to interior  */
 19:   /*       nodes) and                                              */
 20:   /* * N - stands for Neumann (by extension, refers to all local   */
 21:   /*       nodes, interior plus interface).                        */
 22:   /* In some cases, I or D would apply equaly well (e.g. vec1_D).  */

 24:   PetscInt n;                /* number of nodes (interior+interface) in this subdomain */
 25:   PetscInt n_B;              /* number of interface nodes in this subdomain */
 26:   IS       is_B_local,       /* local (sequential) index sets for interface (B) and interior (I) nodes */
 27:            is_I_local,
 28:            is_B_global,
 29:            is_I_global;

 31:   Mat A_II, A_IB,            /* local (sequential) submatrices */
 32:       A_BI, A_BB;
 33:   Mat pA_II;
 34:   Vec D;                     /* diagonal scaling "matrix" (stored as a vector, since it's diagonal) */
 35:   KSP ksp_N,                /* linear solver contexts */
 36:       ksp_D;
 37:   Vec vec1_N,                /* local (sequential) work vectors */
 38:       vec2_N,
 39:       vec1_D,
 40:       vec2_D,
 41:       vec3_D,
 42:       vec4_D,
 43:       vec1_B,
 44:       vec2_B,
 45:       vec3_B,
 46:       vec1_global;

 48:   PetscScalar * work_N;
 49:   VecScatter  N_to_D;             /* scattering context from all local nodes to local interior nodes */
 50:   VecScatter  global_to_D;        /* scattering context from global to local interior nodes */
 51:   VecScatter  N_to_B;             /* scattering context from all local nodes to local interface nodes */
 52:   VecScatter  global_to_B;        /* scattering context from global to local interface nodes */
 53:   PetscBool   pure_neumann;
 54:   PetscScalar scaling_factor;
 55:   PetscBool   use_stiffness_scaling;

 57:   ISLocalToGlobalMapping mapping;
 58:   PetscInt  n_neigh;     /* number of neighbours this subdomain has (INCLUDING the subdomain itself).       */
 59:   PetscInt *neigh;       /* list of neighbouring subdomains                                                 */
 60:   PetscInt *n_shared;    /* n_shared[j] is the number of nodes shared with subdomain neigh[j]               */
 61:   PetscInt **shared;     /* shared[j][i] is the local index of the i-th node shared with subdomain neigh[j] */
 62:   /* It is necessary some consistency in the                                                  */
 63:   /* numbering of the shared edges from each side.                                            */
 64:   /* For instance:                                                                            */
 65:   /*                                                                                          */
 66:   /* +-------+-------+                                                                        */
 67:   /* |   k   |   l   | subdomains k and l are neighbours                                      */
 68:   /* +-------+-------+                                                                        */
 69:   /*                                                                                          */
 70:   /* Let i and j be s.t. proc[k].neigh[i]==l and                                              */
 71:   /*                     proc[l].neigh[j]==k.                                                 */
 72:   /*                                                                                          */
 73:   /* We need:                                                                                 */
 74:   /* proc[k].loc_to_glob(proc[k].shared[i][m]) == proc[l].loc_to_glob(proc[l].shared[j][m])   */
 75:   /* for all 0 <= m < proc[k].n_shared[i], or equiv'ly, for all 0 <= m < proc[l].n_shared[j]  */
 76:   ISLocalToGlobalMapping BtoNmap;
 77:   PetscBool reusesubmatrices;
 78: } PC_IS;

 80: PETSC_EXTERN PetscErrorCode PCISSetUp(PC,PetscBool,PetscBool);
 81: PETSC_EXTERN PetscErrorCode PCISDestroy(PC);
 82: PETSC_EXTERN PetscErrorCode PCISCreate(PC);
 83: PETSC_EXTERN PetscErrorCode PCISApplySchur(PC,Vec,Vec,Vec,Vec,Vec);
 84: PETSC_EXTERN PetscErrorCode PCISScatterArrayNToVecB(PetscScalar*,Vec,InsertMode,ScatterMode,PC);
 85: PETSC_EXTERN PetscErrorCode PCISApplyInvSchur(PC,Vec,Vec,Vec,Vec);

 87: #endif /* __pcis_h */