Actual source code: pcisimpl.h

  1: #pragma once

  3: #include <petsc/private/pcimpl.h>
  4: #include <petsc/private/matisimpl.h>
  5: #include <petscksp.h>

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

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

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

 26:   Mat A_II, A_IB, /* local (sequential) submatrices */
 27:     A_BI, A_BB;
 28:   Mat pA_II;
 29:   Vec D;     /* diagonal scaling "matrix" (stored as a vector, since it's diagonal) */
 30:   KSP ksp_N, /* linear solver contexts */
 31:     ksp_D;
 32:   Vec vec1_N, /* local (sequential) work vectors */
 33:     vec2_N, vec1_D, vec2_D, vec3_D, vec4_D, vec1_B, vec2_B, vec3_B, vec1_global;

 35:   PetscScalar *work_N;
 36:   VecScatter   N_to_D;      /* scattering context from all local nodes to local interior nodes */
 37:   VecScatter   global_to_D; /* scattering context from global to local interior nodes */
 38:   VecScatter   N_to_B;      /* scattering context from all local nodes to local interface nodes */
 39:   VecScatter   global_to_B; /* scattering context from global to local interface nodes */
 40:   PetscBool    pure_neumann;
 41:   PetscScalar  scaling_factor;
 42:   PetscBool    use_stiffness_scaling;

 44:   ISLocalToGlobalMapping mapping;
 45:   PetscInt               n_neigh;  /* number of neighbours this subdomain has (INCLUDING the subdomain itself).       */
 46:   PetscInt              *neigh;    /* list of neighbouring subdomains                                                 */
 47:   PetscInt              *n_shared; /* n_shared[j] is the number of nodes shared with subdomain neigh[j]               */
 48:   PetscInt             **shared;   /* shared[j][i] is the local index of the i-th node shared with subdomain neigh[j] */
 49:   /* It is necessary some consistency in the                                                  */
 50:   /* numbering of the shared edges from each side.                                            */
 51:   /* For instance:                                                                            */
 52:   /*                                                                                          */
 53:   /* +-------+-------+                                                                        */
 54:   /* |   k   |   l   | subdomains k and l are neighbours                                      */
 55:   /* +-------+-------+                                                                        */
 56:   /*                                                                                          */
 57:   /* Let i and j be s.t. proc[k].neigh[i]==l and                                              */
 58:   /*                     proc[l].neigh[j]==k.                                                 */
 59:   /*                                                                                          */
 60:   /* We need:                                                                                 */
 61:   /* proc[k].loc_to_glob(proc[k].shared[i][m]) == proc[l].loc_to_glob(proc[l].shared[j][m])   */
 62:   /* for all 0 <= m < proc[k].n_shared[i], or equiv'ly, for all 0 <= m < proc[l].n_shared[j]  */
 63:   ISLocalToGlobalMapping BtoNmap;
 64:   PetscBool              reusesubmatrices;
 65: } PC_IS;