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