Actual source code: bddcstructs.h

petsc-3.10.5 2019-03-28
Report Typos and Errors

  4:  #include <petscksp.h>
  5:  #include <petscbt.h>

  7: /* special marks for interface graph: they cannot be enums, since special marks should in principle range from -4 to -max_int */
  8: #define PCBDDCGRAPH_NEUMANN_MARK -1
  9: #define PCBDDCGRAPH_DIRICHLET_MARK -2
 10: #define PCBDDCGRAPH_LOCAL_PERIODIC_MARK -3
 11: #define PCBDDCGRAPH_SPECIAL_MARK -4

 13: /* Structure for local graph partitioning */
 14: struct _PCBDDCGraph {
 15:   PetscBool              setupcalled;
 16:   /* graph information */
 17:   ISLocalToGlobalMapping l2gmap;
 18:   PetscInt               nvtxs;
 19:   PetscInt               nvtxs_global;
 20:   PetscBT                touched;
 21:   PetscInt               *count;
 22:   PetscInt               **neighbours_set;
 23:   PetscInt               *subset;
 24:   PetscInt               *which_dof;
 25:   PetscInt               *special_dof;
 26:   PetscInt               custom_minimal_size;
 27:   PetscBool              twodim;
 28:   PetscBool              twodimset;
 29:   PetscBool              has_dirichlet;
 30:   IS                     dirdofs;
 31:   IS                     dirdofsB;
 32:   PetscInt               commsizelimit;
 33:   PetscInt               maxcount;
 34:   /* data for connected components */
 35:   PetscInt               ncc;
 36:   PetscInt               *cptr;
 37:   PetscInt               *queue;
 38:   PetscBool              queue_sorted;
 39:   /* data for interface subsets */
 40:   PetscInt               n_subsets;
 41:   PetscInt               *subset_size;
 42:   PetscInt               **subset_idxs;
 43:   PetscInt               *subset_ncc;
 44:   PetscInt               *subset_ref_node;
 45:   /* data for periodic dofs */
 46:   PetscInt               *mirrors;
 47:   PetscInt               **mirrors_set;
 48:   /* placeholders for connectivity relation between dofs */
 49:   PetscInt               nvtxs_csr;
 50:   PetscInt               *xadj;
 51:   PetscInt               *adjncy;
 52:   PetscBool              freecsr;
 53:   /* data for local subdomains (if any have been detected)
 54:      these are not intended to be exposed */
 55:   PetscInt               n_local_subs;
 56:   PetscInt               *local_subs;
 57:   /* coordinates (for corner detection) */
 58:   PetscBool              active_coords;
 59:   PetscBool              cloc;
 60:   PetscInt               cdim,cnloc;
 61:   PetscReal*             coords;

 63: };
 64: typedef struct _PCBDDCGraph *PCBDDCGraph;

 66: /* Wrap to MatFactor solver in Schur complement mode. Provides
 67:    - standalone solver for interior variables
 68:    - forward and backward substitutions for correction solver
 69: */
 70: /* It assumes that interior variables are a contiguous set starting from 0 */
 71: struct _PCBDDCReuseSolvers {
 72:   /* the factored matrix obtained from MatGetFactor(...,solver_package,...) */
 73:   Mat        F;
 74:   /* placeholders for the solution and rhs on the whole set of dofs of A (size local_dofs - local_vertices)*/
 75:   Vec        sol;
 76:   Vec        rhs;
 77:   /* */
 78:   PetscBool  has_vertices;
 79:   /* shell PCs to handle interior/correction solvers */
 80:   PC         interior_solver;
 81:   PC         correction_solver;
 82:   IS         is_R;
 83:   /* objects to handle Schur complement solution */
 84:   Vec        rhs_B;
 85:   Vec        sol_B;
 86:   IS         is_B;
 87:   VecScatter correction_scatter_B;
 88:   /* handle benign trick without change of basis on pressures */
 89:   PetscInt    benign_n;
 90:   IS          *benign_zerodiag_subs;
 91:   PetscScalar *benign_save_vals;
 92:   Mat         benign_csAIB;
 93:   Mat         benign_AIIm1ones;
 94:   Vec         benign_corr_work;
 95:   Vec         benign_dummy_schur_vec;
 96: };
 97: typedef struct _PCBDDCReuseSolvers *PCBDDCReuseSolvers;

 99: /* structure to handle Schur complements on subsets */
100: struct _PCBDDCSubSchurs {
101:   /* local Neumann matrix */
102:   Mat A;
103:   /* local Schur complement */
104:   Mat S;
105:   /* index sets */
106:   IS  is_I;
107:   IS  is_B;
108:   /* whether Schur complements are explicitly computed with or not */
109:   char      mat_solver_type[64];
110:   PetscBool schur_explicit;
111:   /* matrices cointained explicit schur complements cat together */
112:   /* note that AIJ format is used but the values are inserted as in column major ordering */
113:   Mat S_Ej_all;
114:   Mat sum_S_Ej_all;
115:   Mat sum_S_Ej_inv_all;
116:   Mat sum_S_Ej_tilda_all;
117:   IS  is_Ej_all;
118:   IS  is_vertices;
119:   IS  is_dir;
120:   /* l2g maps */
121:   ISLocalToGlobalMapping l2gmap;
122:   ISLocalToGlobalMapping BtoNmap;
123:   /* number of local subproblems */
124:   PetscInt n_subs;
125:   /* connected components */
126:   IS*      is_subs;
127:   PetscBT  is_edge;
128:   /* mat flags */
129:   PetscBool is_symmetric;
130:   PetscBool is_hermitian;
131:   PetscBool is_posdef;
132:   /* data structure to reuse MatFactor with Schur solver */
133:   PCBDDCReuseSolvers reuse_solver;
134:   /* change of variables */
135:   KSP       *change;
136:   IS        *change_primal_sub;
137:   PetscBool change_with_qr;
138:   /* prefix */
139:   char      *prefix;
140:   /* */
141:   PetscBool restrict_comm;
142:   /* debug */
143:   PetscBool debug;
144: };
145: typedef struct _PCBDDCSubSchurs *PCBDDCSubSchurs;

147: /* Structure for deluxe scaling */
148: struct _PCBDDCDeluxeScaling {
149:   /* simple scaling on selected dofs (i.e. primal vertices and nodes on interface dirichlet boundaries) */
150:   PetscInt        n_simple;
151:   PetscInt*       idx_simple_B;
152:   /* handle deluxe problems  */
153:   PetscInt        seq_n;
154:   PetscScalar     *workspace;
155:   VecScatter      *seq_scctx;
156:   Vec             *seq_work1;
157:   Vec             *seq_work2;
158:   Mat             *seq_mat;
159:   Mat             *seq_mat_inv_sum;
160:   KSP             *change;
161:   PetscBool       change_with_qr;
162: };
163: typedef struct _PCBDDCDeluxeScaling *PCBDDCDeluxeScaling;

165: /* inexact solvers with nullspace correction */
166: struct _NullSpaceCorrection_ctx {
167:   Mat         basis_mat;
168:   Mat         Kbasis_mat;
169:   Mat         Lbasis_mat;
170:   PC          local_pc;
171:   Vec         work_small_1;
172:   Vec         work_small_2;
173:   Vec         work_full_1;
174:   Vec         work_full_2;
175:   PetscBool   apply_scaling;
176:   PetscScalar scale;
177: };
178: typedef struct _NullSpaceCorrection_ctx *NullSpaceCorrection_ctx;

180: /* MatShell context for benign mat mults */
181: struct _PCBDDCBenignMatMult_ctx {
182:   Mat         A;
183:   PetscInt    benign_n;
184:   IS          *benign_zerodiag_subs;
185:   PetscScalar *work;
186:   PetscBool   apply_left;
187:   PetscBool   apply_right;
188:   PetscBool   apply_p0;
189:   PetscBool   free;
190: };
191: typedef struct _PCBDDCBenignMatMult_ctx *PCBDDCBenignMatMult_ctx;

193: /* feti-dp mat */
194: struct _FETIDPMat_ctx {
195:   PetscInt   n;               /* local number of rows */
196:   PetscInt   N;               /* global number of rows */
197:   PetscInt   n_lambda;        /* global number of multipliers */
198:   Vec        lambda_local;
199:   Vec        temp_solution_B;
200:   Vec        temp_solution_D;
201:   Mat        B_delta;
202:   Mat        B_Ddelta;
203:   PetscBool  deluxe_nonred;
204:   VecScatter l2g_lambda;
205:   PC         pc;
206:   PetscBool  fully_redundant;
207:   /* saddle point */
208:   VecScatter l2g_lambda_only;
209:   Mat        B_BB;
210:   Mat        B_BI;
211:   Mat        Bt_BB;
212:   Mat        Bt_BI;
213:   Mat        C;
214:   VecScatter l2g_p;
215:   VecScatter g2g_p;
216:   Vec        vP;
217:   Vec        xPg;
218:   Vec        yPg;
219:   Vec        rhs_flip;
220:   IS         pressure;
221:   IS         lagrange;
222: };
223: typedef struct _FETIDPMat_ctx *FETIDPMat_ctx;

225: /* feti-dp preconditioner */
226: struct _FETIDPPC_ctx {
227:   Mat        S_j;
228:   Vec        lambda_local;
229:   Mat        B_Ddelta;
230:   VecScatter l2g_lambda;
231:   PC         pc;
232:   /* saddle point */
233:   Vec        xPg;
234:   Vec        yPg;
235: };
236: typedef struct _FETIDPPC_ctx *FETIDPPC_ctx;

238: struct _BDdelta_DN {
239:   Mat BD;
240:   KSP kBD;
241:   Vec work;
242: };
243: typedef struct _BDdelta_DN *BDdelta_DN;

245: /* Schur interface preconditioner */
246: struct _BDDCIPC_ctx {
247:   VecScatter g2l;
248:   PC         bddc;
249: };
250: typedef struct _BDDCIPC_ctx *BDDCIPC_ctx;

252: #endif