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