xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision daa037dfd3c3bec8dc8659548d2b20b07c1dc6de)
1 
2 #include <../src/ksp/pc/impls/is/pcis.h> /*I "petscpc.h" I*/
3 
4 static PetscErrorCode PCISSetUseStiffnessScaling_IS(PC pc, PetscBool use)
5 {
6   PC_IS *pcis = (PC_IS*)pc->data;
7 
8   PetscFunctionBegin;
9   pcis->use_stiffness_scaling = use;
10   PetscFunctionReturn(0);
11 }
12 
13 /*@
14  PCISSetUseStiffnessScaling - Tells PCIS to construct partition of unity using
15                               local matrices' diagonal.
16 
17    Not collective
18 
19    Input Parameters:
20 +  pc - the preconditioning context
21 -  use - whether or not pcis use matrix diagonal to build partition of unity.
22 
23    Level: intermediate
24 
25    Notes:
26 
27 .seealso: PCBDDC
28 @*/
29 PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use)
30 {
31   PetscFunctionBegin;
32   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
33   PetscValidLogicalCollectiveBool(pc,use,2);
34   PetscTryMethod(pc,"PCISSetUseStiffnessScaling_C",(PC,PetscBool),(pc,use));
35   PetscFunctionReturn(0);
36 }
37 
38 static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors)
39 {
40   PC_IS          *pcis = (PC_IS*)pc->data;
41 
42   PetscFunctionBegin;
43   PetscCall(PetscObjectReference((PetscObject)scaling_factors));
44   PetscCall(VecDestroy(&pcis->D));
45   pcis->D = scaling_factors;
46   if (pc->setupcalled) {
47     PetscInt sn;
48 
49     PetscCall(VecGetSize(pcis->D,&sn));
50     if (sn == pcis->n) {
51       PetscCall(VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
52       PetscCall(VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
53       PetscCall(VecDestroy(&pcis->D));
54       PetscCall(VecDuplicate(pcis->vec1_B,&pcis->D));
55       PetscCall(VecCopy(pcis->vec1_B,pcis->D));
56     } else PetscCheckFalse(sn != pcis->n_B,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Invalid size for scaling vector. Expected %D (or full %D), found %D",pcis->n_B,pcis->n,sn);
57   }
58   PetscFunctionReturn(0);
59 }
60 
61 /*@
62  PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS.
63 
64    Not collective
65 
66    Input Parameters:
67 +  pc - the preconditioning context
68 -  scaling_factors - scaling factors for the subdomain
69 
70    Level: intermediate
71 
72    Notes:
73    Intended to use with jumping coefficients cases.
74 
75 .seealso: PCBDDC
76 @*/
77 PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors)
78 {
79   PetscFunctionBegin;
80   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
81   PetscValidHeaderSpecific(scaling_factors,VEC_CLASSID,2);
82   PetscTryMethod(pc,"PCISSetSubdomainDiagonalScaling_C",(PC,Vec),(pc,scaling_factors));
83   PetscFunctionReturn(0);
84 }
85 
86 static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal)
87 {
88   PC_IS *pcis = (PC_IS*)pc->data;
89 
90   PetscFunctionBegin;
91   pcis->scaling_factor = scal;
92   if (pcis->D) {
93 
94     PetscCall(VecSet(pcis->D,pcis->scaling_factor));
95   }
96   PetscFunctionReturn(0);
97 }
98 
99 /*@
100  PCISSetSubdomainScalingFactor - Set scaling factor for PCIS.
101 
102    Not collective
103 
104    Input Parameters:
105 +  pc - the preconditioning context
106 -  scal - scaling factor for the subdomain
107 
108    Level: intermediate
109 
110    Notes:
111    Intended to use with jumping coefficients cases.
112 
113 .seealso: PCBDDC
114 @*/
115 PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal)
116 {
117   PetscFunctionBegin;
118   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
119   PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));
120   PetscFunctionReturn(0);
121 }
122 
123 /* -------------------------------------------------------------------------- */
124 /*
125    PCISSetUp -
126 */
127 PetscErrorCode  PCISSetUp(PC pc, PetscBool computematrices, PetscBool computesolvers)
128 {
129   PC_IS          *pcis  = (PC_IS*)(pc->data);
130   Mat_IS         *matis;
131   MatReuse       reuse;
132   PetscErrorCode ierr;
133   PetscBool      flg,issbaij;
134 
135   PetscFunctionBegin;
136   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&flg));
137   PetscCheck(flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Requires preconditioning matrix of type MATIS");
138   matis = (Mat_IS*)pc->pmat->data;
139   if (pc->useAmat) {
140     PetscCall(PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&flg));
141     PetscCheck(flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Requires linear system matrix of type MATIS");
142   }
143 
144   /* first time creation, get info on substructuring */
145   if (!pc->setupcalled) {
146     PetscInt    n_I;
147     PetscInt    *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global;
148     PetscBT     bt;
149     PetscInt    i,j;
150 
151     /* get info on mapping */
152     PetscCall(PetscObjectReference((PetscObject)matis->rmapping));
153     PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
154     pcis->mapping = matis->rmapping;
155     PetscCall(ISLocalToGlobalMappingGetSize(pcis->mapping,&pcis->n));
156     PetscCall(ISLocalToGlobalMappingGetInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared)));
157 
158     /* Identifying interior and interface nodes, in local numbering */
159     PetscCall(PetscBTCreate(pcis->n,&bt));
160     for (i=0;i<pcis->n_neigh;i++)
161       for (j=0;j<pcis->n_shared[i];j++) {
162         PetscCall(PetscBTSet(bt,pcis->shared[i][j]));
163       }
164 
165     /* Creating local and global index sets for interior and inteface nodes. */
166     PetscCall(PetscMalloc1(pcis->n,&idx_I_local));
167     PetscCall(PetscMalloc1(pcis->n,&idx_B_local));
168     for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) {
169       if (!PetscBTLookup(bt,i)) {
170         idx_I_local[n_I] = i;
171         n_I++;
172       } else {
173         idx_B_local[pcis->n_B] = i;
174         pcis->n_B++;
175       }
176     }
177 
178     /* Getting the global numbering */
179     idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
180     idx_I_global = idx_B_local + pcis->n_B;
181     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping,pcis->n_B,idx_B_local,idx_B_global));
182     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping,n_I,idx_I_local,idx_I_global));
183 
184     /* Creating the index sets */
185     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local));
186     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc),pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global));
187     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local));
188     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc),n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global));
189 
190     /* Freeing memory */
191     PetscCall(PetscFree(idx_B_local));
192     PetscCall(PetscFree(idx_I_local));
193     PetscCall(PetscBTDestroy(&bt));
194 
195     /* Creating work vectors and arrays */
196     PetscCall(VecDuplicate(matis->x,&pcis->vec1_N));
197     PetscCall(VecDuplicate(pcis->vec1_N,&pcis->vec2_N));
198     PetscCall(VecCreate(PETSC_COMM_SELF,&pcis->vec1_D));
199     PetscCall(VecSetSizes(pcis->vec1_D,pcis->n-pcis->n_B,PETSC_DECIDE));
200     PetscCall(VecSetType(pcis->vec1_D,((PetscObject)pcis->vec1_N)->type_name));
201     PetscCall(VecDuplicate(pcis->vec1_D,&pcis->vec2_D));
202     PetscCall(VecDuplicate(pcis->vec1_D,&pcis->vec3_D));
203     PetscCall(VecDuplicate(pcis->vec1_D,&pcis->vec4_D));
204     PetscCall(VecCreate(PETSC_COMM_SELF,&pcis->vec1_B));
205     PetscCall(VecSetSizes(pcis->vec1_B,pcis->n_B,PETSC_DECIDE));
206     PetscCall(VecSetType(pcis->vec1_B,((PetscObject)pcis->vec1_N)->type_name));
207     PetscCall(VecDuplicate(pcis->vec1_B,&pcis->vec2_B));
208     PetscCall(VecDuplicate(pcis->vec1_B,&pcis->vec3_B));
209     PetscCall(MatCreateVecs(pc->pmat,&pcis->vec1_global,NULL));
210     PetscCall(PetscMalloc1(pcis->n,&pcis->work_N));
211     /* scaling vector */
212     if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */
213       PetscCall(VecDuplicate(pcis->vec1_B,&pcis->D));
214       PetscCall(VecSet(pcis->D,pcis->scaling_factor));
215     }
216 
217     /* Creating the scatter contexts */
218     PetscCall(VecScatterCreate(pcis->vec1_N,pcis->is_I_local,pcis->vec1_D,(IS)0,&pcis->N_to_D));
219     PetscCall(VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D));
220     PetscCall(VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B));
221     PetscCall(VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B));
222 
223     /* map from boundary to local */
224     PetscCall(ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcis->BtoNmap));
225   }
226 
227   {
228     PetscInt sn;
229 
230     PetscCall(VecGetSize(pcis->D,&sn));
231     if (sn == pcis->n) {
232       PetscCall(VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
233       PetscCall(VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
234       PetscCall(VecDestroy(&pcis->D));
235       PetscCall(VecDuplicate(pcis->vec1_B,&pcis->D));
236       PetscCall(VecCopy(pcis->vec1_B,pcis->D));
237     } else PetscCheckFalse(sn != pcis->n_B,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Invalid size for scaling vector. Expected %D (or full %D), found %D",pcis->n_B,pcis->n,sn);
238   }
239 
240   /*
241     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
242     is such that interior nodes come first than the interface ones, we have
243 
244         [ A_II | A_IB ]
245     A = [------+------]
246         [ A_BI | A_BB ]
247   */
248   if (computematrices) {
249     PetscBool amat = (PetscBool)(pc->mat != pc->pmat && pc->useAmat);
250     PetscInt  bs,ibs;
251 
252     reuse = MAT_INITIAL_MATRIX;
253     if (pcis->reusesubmatrices && pc->setupcalled) {
254       if (pc->flag == SAME_NONZERO_PATTERN) {
255         reuse = MAT_REUSE_MATRIX;
256       } else {
257         reuse = MAT_INITIAL_MATRIX;
258       }
259     }
260     if (reuse == MAT_INITIAL_MATRIX) {
261       PetscCall(MatDestroy(&pcis->A_II));
262       PetscCall(MatDestroy(&pcis->pA_II));
263       PetscCall(MatDestroy(&pcis->A_IB));
264       PetscCall(MatDestroy(&pcis->A_BI));
265       PetscCall(MatDestroy(&pcis->A_BB));
266     }
267 
268     PetscCall(ISLocalToGlobalMappingGetBlockSize(pcis->mapping,&ibs));
269     PetscCall(MatGetBlockSize(matis->A,&bs));
270     PetscCall(MatCreateSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->pA_II));
271     if (amat) {
272       Mat_IS *amatis = (Mat_IS*)pc->mat->data;
273       PetscCall(MatCreateSubMatrix(amatis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II));
274     } else {
275       PetscCall(PetscObjectReference((PetscObject)pcis->pA_II));
276       PetscCall(MatDestroy(&pcis->A_II));
277       pcis->A_II = pcis->pA_II;
278     }
279     PetscCall(MatSetBlockSize(pcis->A_II,bs == ibs ? bs : 1));
280     PetscCall(MatSetBlockSize(pcis->pA_II,bs == ibs ? bs : 1));
281     PetscCall(MatCreateSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB));
282     PetscCall(PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij));
283     if (!issbaij) {
284       PetscCall(MatCreateSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB));
285       PetscCall(MatCreateSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI));
286     } else {
287       Mat newmat;
288 
289       PetscCall(MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&newmat));
290       PetscCall(MatCreateSubMatrix(newmat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB));
291       PetscCall(MatCreateSubMatrix(newmat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI));
292       PetscCall(MatDestroy(&newmat));
293     }
294     PetscCall(MatSetBlockSize(pcis->A_BB,bs == ibs ? bs : 1));
295   }
296 
297   /* Creating scaling vector D */
298   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL));
299   if (pcis->use_stiffness_scaling) {
300     PetscScalar *a;
301     PetscInt    i,n;
302 
303     if (pcis->A_BB) {
304       PetscCall(MatGetDiagonal(pcis->A_BB,pcis->D));
305     } else {
306       PetscCall(MatGetDiagonal(matis->A,pcis->vec1_N));
307       PetscCall(VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD));
308       PetscCall(VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD));
309     }
310     PetscCall(VecAbs(pcis->D));
311     PetscCall(VecGetLocalSize(pcis->D,&n));
312     PetscCall(VecGetArray(pcis->D,&a));
313     for (i=0;i<n;i++) if (PetscAbsScalar(a[i])<PETSC_SMALL) a[i] = 1.0;
314     PetscCall(VecRestoreArray(pcis->D,&a));
315   }
316   PetscCall(VecSet(pcis->vec1_global,0.0));
317   PetscCall(VecScatterBegin(pcis->global_to_B,pcis->D,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
318   PetscCall(VecScatterEnd(pcis->global_to_B,pcis->D,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
319   PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
320   PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
321   PetscCall(VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B));
322   /* See historical note 01, at the bottom of this file. */
323 
324   /* Creating the KSP contexts for the local Dirichlet and Neumann problems */
325   if (computesolvers) {
326     PC pc_ctx;
327 
328     pcis->pure_neumann = matis->pure_neumann;
329     /* Dirichlet */
330     PetscCall(KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D));
331     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_D,pc->erroriffailure));
332     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1));
333     PetscCall(KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II));
334     PetscCall(KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_"));
335     PetscCall(KSPGetPC(pcis->ksp_D,&pc_ctx));
336     PetscCall(PCSetType(pc_ctx,PCLU));
337     PetscCall(KSPSetType(pcis->ksp_D,KSPPREONLY));
338     PetscCall(KSPSetFromOptions(pcis->ksp_D));
339     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
340     PetscCall(KSPSetUp(pcis->ksp_D));
341     /* Neumann */
342     PetscCall(KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N));
343     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_N,pc->erroriffailure));
344     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1));
345     PetscCall(KSPSetOperators(pcis->ksp_N,matis->A,matis->A));
346     PetscCall(KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_"));
347     PetscCall(KSPGetPC(pcis->ksp_N,&pc_ctx));
348     PetscCall(PCSetType(pc_ctx,PCLU));
349     PetscCall(KSPSetType(pcis->ksp_N,KSPPREONLY));
350     PetscCall(KSPSetFromOptions(pcis->ksp_N));
351     {
352       PetscBool damp_fixed                    = PETSC_FALSE,
353                 remove_nullspace_fixed        = PETSC_FALSE,
354                 set_damping_factor_floating   = PETSC_FALSE,
355                 not_damp_floating             = PETSC_FALSE,
356                 not_remove_nullspace_floating = PETSC_FALSE;
357       PetscReal fixed_factor,
358                 floating_factor;
359 
360       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed));
361       if (!damp_fixed) fixed_factor = 0.0;
362       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL));
363 
364       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL));
365 
366       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",
367                               &floating_factor,&set_damping_factor_floating);PetscCall(ierr);
368       if (!set_damping_factor_floating) floating_factor = 0.0;
369       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL));
370       if (!set_damping_factor_floating) floating_factor = 1.e-12;
371 
372       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",&not_damp_floating,NULL));
373 
374       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",&not_remove_nullspace_floating,NULL));
375 
376       if (pcis->pure_neumann) {  /* floating subdomain */
377         if (!(not_damp_floating)) {
378           PetscCall(PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO));
379           PetscCall(PCFactorSetShiftAmount(pc_ctx,floating_factor));
380         }
381         if (!(not_remove_nullspace_floating)) {
382           MatNullSpace nullsp;
383           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp));
384           PetscCall(MatSetNullSpace(matis->A,nullsp));
385           PetscCall(MatNullSpaceDestroy(&nullsp));
386         }
387       } else {  /* fixed subdomain */
388         if (damp_fixed) {
389           PetscCall(PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO));
390           PetscCall(PCFactorSetShiftAmount(pc_ctx,floating_factor));
391         }
392         if (remove_nullspace_fixed) {
393           MatNullSpace nullsp;
394           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp));
395           PetscCall(MatSetNullSpace(matis->A,nullsp));
396           PetscCall(MatNullSpaceDestroy(&nullsp));
397         }
398       }
399     }
400     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
401     PetscCall(KSPSetUp(pcis->ksp_N));
402   }
403   PetscFunctionReturn(0);
404 }
405 
406 /* -------------------------------------------------------------------------- */
407 /*
408    PCISDestroy -
409 */
410 PetscErrorCode  PCISDestroy(PC pc)
411 {
412   PC_IS          *pcis = (PC_IS*)(pc->data);
413 
414   PetscFunctionBegin;
415   PetscCall(ISDestroy(&pcis->is_B_local));
416   PetscCall(ISDestroy(&pcis->is_I_local));
417   PetscCall(ISDestroy(&pcis->is_B_global));
418   PetscCall(ISDestroy(&pcis->is_I_global));
419   PetscCall(MatDestroy(&pcis->A_II));
420   PetscCall(MatDestroy(&pcis->pA_II));
421   PetscCall(MatDestroy(&pcis->A_IB));
422   PetscCall(MatDestroy(&pcis->A_BI));
423   PetscCall(MatDestroy(&pcis->A_BB));
424   PetscCall(VecDestroy(&pcis->D));
425   PetscCall(KSPDestroy(&pcis->ksp_N));
426   PetscCall(KSPDestroy(&pcis->ksp_D));
427   PetscCall(VecDestroy(&pcis->vec1_N));
428   PetscCall(VecDestroy(&pcis->vec2_N));
429   PetscCall(VecDestroy(&pcis->vec1_D));
430   PetscCall(VecDestroy(&pcis->vec2_D));
431   PetscCall(VecDestroy(&pcis->vec3_D));
432   PetscCall(VecDestroy(&pcis->vec4_D));
433   PetscCall(VecDestroy(&pcis->vec1_B));
434   PetscCall(VecDestroy(&pcis->vec2_B));
435   PetscCall(VecDestroy(&pcis->vec3_B));
436   PetscCall(VecDestroy(&pcis->vec1_global));
437   PetscCall(VecScatterDestroy(&pcis->global_to_D));
438   PetscCall(VecScatterDestroy(&pcis->N_to_B));
439   PetscCall(VecScatterDestroy(&pcis->N_to_D));
440   PetscCall(VecScatterDestroy(&pcis->global_to_B));
441   PetscCall(PetscFree(pcis->work_N));
442   if (pcis->n_neigh > -1) {
443     PetscCall(ISLocalToGlobalMappingRestoreInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared)));
444   }
445   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
446   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->BtoNmap));
447   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL));
448   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL));
449   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL));
450   PetscFunctionReturn(0);
451 }
452 
453 /* -------------------------------------------------------------------------- */
454 /*
455    PCISCreate -
456 */
457 PetscErrorCode  PCISCreate(PC pc)
458 {
459   PC_IS          *pcis = (PC_IS*)(pc->data);
460 
461   PetscFunctionBegin;
462   pcis->n_neigh          = -1;
463   pcis->scaling_factor   = 1.0;
464   pcis->reusesubmatrices = PETSC_TRUE;
465   /* composing functions */
466   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS));
467   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS));
468   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS));
469   PetscFunctionReturn(0);
470 }
471 
472 /* -------------------------------------------------------------------------- */
473 /*
474    PCISApplySchur -
475 
476    Input parameters:
477 .  pc - preconditioner context
478 .  v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
479 
480    Output parameters:
481 .  vec1_B - result of Schur complement applied to chunk
482 .  vec2_B - garbage (used as work space), or null (and v is used as workspace)
483 .  vec1_D - garbage (used as work space)
484 .  vec2_D - garbage (used as work space)
485 
486 */
487 PetscErrorCode  PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
488 {
489   PC_IS          *pcis = (PC_IS*)(pc->data);
490 
491   PetscFunctionBegin;
492   if (!vec2_B) vec2_B = v;
493 
494   PetscCall(MatMult(pcis->A_BB,v,vec1_B));
495   PetscCall(MatMult(pcis->A_IB,v,vec1_D));
496   PetscCall(KSPSolve(pcis->ksp_D,vec1_D,vec2_D));
497   PetscCall(KSPCheckSolve(pcis->ksp_D,pc,vec2_D));
498   PetscCall(MatMult(pcis->A_BI,vec2_D,vec2_B));
499   PetscCall(VecAXPY(vec1_B,-1.0,vec2_B));
500   PetscFunctionReturn(0);
501 }
502 
503 /* -------------------------------------------------------------------------- */
504 /*
505    PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
506    including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
507    mode.
508 
509    Input parameters:
510 .  pc - preconditioner context
511 .  array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
512 .  v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
513 
514    Output parameter:
515 .  array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
516 .  v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
517 
518    Notes:
519    The entries in the array that do not correspond to interface nodes remain unaltered.
520 */
521 PetscErrorCode  PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc)
522 {
523   PetscInt       i;
524   const PetscInt *idex;
525   PetscScalar    *array_B;
526   PC_IS          *pcis = (PC_IS*)(pc->data);
527 
528   PetscFunctionBegin;
529   PetscCall(VecGetArray(v_B,&array_B));
530   PetscCall(ISGetIndices(pcis->is_B_local,&idex));
531 
532   if (smode == SCATTER_FORWARD) {
533     if (imode == INSERT_VALUES) {
534       for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]];
535     } else {  /* ADD_VALUES */
536       for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]];
537     }
538   } else {  /* SCATTER_REVERSE */
539     if (imode == INSERT_VALUES) {
540       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i];
541     } else {  /* ADD_VALUES */
542       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i];
543     }
544   }
545   PetscCall(ISRestoreIndices(pcis->is_B_local,&idex));
546   PetscCall(VecRestoreArray(v_B,&array_B));
547   PetscFunctionReturn(0);
548 }
549 
550 /* -------------------------------------------------------------------------- */
551 /*
552    PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
553    More precisely, solves the problem:
554                                         [ A_II  A_IB ] [ . ]   [ 0 ]
555                                         [            ] [   ] = [   ]
556                                         [ A_BI  A_BB ] [ x ]   [ b ]
557 
558    Input parameters:
559 .  pc - preconditioner context
560 .  b - vector of local interface nodes (including ghosts)
561 
562    Output parameters:
563 .  x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
564        complement to b
565 .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
566 .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
567 
568 */
569 PetscErrorCode  PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
570 {
571   PC_IS          *pcis = (PC_IS*)(pc->data);
572 
573   PetscFunctionBegin;
574   /*
575     Neumann solvers.
576     Applying the inverse of the local Schur complement, i.e, solving a Neumann
577     Problem with zero at the interior nodes of the RHS and extracting the interface
578     part of the solution. inverse Schur complement is applied to b and the result
579     is stored in x.
580   */
581   /* Setting the RHS vec1_N */
582   PetscCall(VecSet(vec1_N,0.0));
583   PetscCall(VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE));
584   PetscCall(VecScatterEnd  (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE));
585   /* Checking for consistency of the RHS */
586   {
587     PetscBool flg = PETSC_FALSE;
588     PetscCall(PetscOptionsGetBool(NULL,NULL,"-pc_is_check_consistency",&flg,NULL));
589     if (flg) {
590       PetscScalar average;
591       PetscViewer viewer;
592       PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer));
593 
594       PetscCall(VecSum(vec1_N,&average));
595       average = average / ((PetscReal)pcis->n);
596       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
597       if (pcis->pure_neumann) {
598         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average)));
599       } else {
600         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed.    Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average)));
601       }
602       PetscCall(PetscViewerFlush(viewer));
603       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
604     }
605   }
606   /* Solving the system for vec2_N */
607   PetscCall(KSPSolve(pcis->ksp_N,vec1_N,vec2_N));
608   PetscCall(KSPCheckSolve(pcis->ksp_N,pc,vec2_N));
609   /* Extracting the local interface vector out of the solution */
610   PetscCall(VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD));
611   PetscCall(VecScatterEnd  (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD));
612   PetscFunctionReturn(0);
613 }
614