xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision 3975b054dc2ded1a903f339d6a1c859863e3c05f)
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, PetscBool computesolvers)
132 {
133   PC_IS          *pcis  = (PC_IS*)(pc->data);
134   Mat_IS         *matis;
135   MatReuse       reuse;
136   PetscErrorCode ierr;
137   PetscBool      flg,issbaij;
138   Vec            counter;
139 
140   PetscFunctionBegin;
141   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&flg);CHKERRQ(ierr);
142   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS");
143   matis = (Mat_IS*)pc->pmat->data;
144 
145   /* first time creation, get info on substructuring */
146   if (!pc->setupcalled) {
147     PetscInt    n_I;
148     PetscInt    *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global;
149     PetscInt    *array;
150     PetscInt    i,j;
151 
152     /* get info on mapping */
153     ierr = PetscObjectReference((PetscObject)matis->mapping);CHKERRQ(ierr);
154     ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr);
155     pcis->mapping = matis->mapping;
156     ierr = ISLocalToGlobalMappingGetSize(pcis->mapping,&pcis->n);CHKERRQ(ierr);
157     ierr = ISLocalToGlobalMappingGetInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
158 
159     /* Identifying interior and interface nodes, in local numbering */
160     ierr = PetscMalloc1(pcis->n,&array);CHKERRQ(ierr);
161     ierr = PetscMemzero(array,pcis->n*sizeof(PetscInt));CHKERRQ(ierr);
162     for (i=0;i<pcis->n_neigh;i++)
163       for (j=0;j<pcis->n_shared[i];j++)
164           array[pcis->shared[i][j]] += 1;
165 
166     /* Creating local and global index sets for interior and inteface nodes. */
167     ierr = PetscMalloc1(pcis->n,&idx_I_local);CHKERRQ(ierr);
168     ierr = PetscMalloc1(pcis->n,&idx_B_local);CHKERRQ(ierr);
169     for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) {
170       if (!array[i]) {
171         idx_I_local[n_I] = i;
172         n_I++;
173       } else {
174         idx_B_local[pcis->n_B] = i;
175         pcis->n_B++;
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     ierr         = ISLocalToGlobalMappingApply(pcis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr);
182     ierr         = ISLocalToGlobalMappingApply(pcis->mapping,n_I,idx_I_local,idx_I_global);CHKERRQ(ierr);
183 
184     /* Creating the index sets */
185     ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr);
186     ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr);
187     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr);
188     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr);
189 
190     /* Freeing memory */
191     ierr = PetscFree(idx_B_local);CHKERRQ(ierr);
192     ierr = PetscFree(idx_I_local);CHKERRQ(ierr);
193     ierr = PetscFree(array);CHKERRQ(ierr);
194 
195     /* Creating work vectors and arrays */
196     ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr);
197     ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr);
198     ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr);
199     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr);
200     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr);
201     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec4_D);CHKERRQ(ierr);
202     ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr);
203     ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr);
204     ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr);
205     ierr = MatCreateVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr);
206     ierr = PetscMalloc1(pcis->n,&pcis->work_N);CHKERRQ(ierr);
207     /* scaling vector */
208     ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr);
209 
210     /* Creating the scatter contexts */
211     ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr);
212     ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr);
213     ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr);
214 
215     /* map from boundary to local */
216     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcis->BtoNmap);CHKERRQ(ierr);
217   }
218 
219   /*
220     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
221     is such that interior nodes come first than the interface ones, we have
222 
223         [ A_II | A_IB ]
224     A = [------+------]
225         [ A_BI | A_BB ]
226   */
227   reuse = MAT_INITIAL_MATRIX;
228   if (pcis->reusesubmatrices && pc->setupcalled) {
229     if (pc->flag == SAME_NONZERO_PATTERN) {
230       reuse = MAT_REUSE_MATRIX;
231     } else {
232       reuse = MAT_INITIAL_MATRIX;
233     }
234   }
235   if (reuse == MAT_INITIAL_MATRIX) {
236     ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
237     ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
238     ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
239     ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
240   }
241 
242   ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr);
243   ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr);
244   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
245   if (!issbaij) {
246     ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
247     ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
248   } else {
249     Mat newmat;
250     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&newmat);CHKERRQ(ierr);
251     ierr = MatGetSubMatrix(newmat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
252     ierr = MatGetSubMatrix(newmat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
253     ierr = MatDestroy(&newmat);CHKERRQ(ierr);
254   }
255 
256   /* Creating scaling "matrix" D */
257   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr);
258   if (!pcis->use_stiffness_scaling) {
259     ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr);
260   } else {
261     ierr = MatGetDiagonal(matis->A,pcis->vec1_N);CHKERRQ(ierr);
262     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
263     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
264   }
265   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
266   ierr = MatCreateVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */
267   ierr = VecSet(counter,0.0);CHKERRQ(ierr);
268   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
269   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
270   ierr = VecScatterBegin(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
271   ierr = VecScatterEnd(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
272   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
273   ierr = VecDestroy(&counter);CHKERRQ(ierr);
274 
275   /* See historical note 01, at the bottom of this file. */
276 
277   /* Creating the KSP contexts for the local Dirichlet and Neumann problems */
278   if (computesolvers) {
279     PC pc_ctx;
280 
281     pcis->pure_neumann = matis->pure_neumann;
282     /* Dirichlet */
283     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr);
284     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr);
285     ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II);CHKERRQ(ierr);
286     ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr);
287     ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr);
288     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
289     ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr);
290     ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr);
291     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
292     ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr);
293     /* Neumann */
294     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr);
295     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr);
296     ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A);CHKERRQ(ierr);
297     ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr);
298     ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr);
299     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
300     ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr);
301     ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr);
302     {
303       PetscBool damp_fixed                    = PETSC_FALSE,
304                 remove_nullspace_fixed        = PETSC_FALSE,
305                 set_damping_factor_floating   = PETSC_FALSE,
306                 not_damp_floating             = PETSC_FALSE,
307                 not_remove_nullspace_floating = PETSC_FALSE;
308       PetscReal fixed_factor,
309                 floating_factor;
310 
311       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr);
312       if (!damp_fixed) fixed_factor = 0.0;
313       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL);CHKERRQ(ierr);
314 
315       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL);CHKERRQ(ierr);
316 
317       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",
318                               &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr);
319       if (!set_damping_factor_floating) floating_factor = 0.0;
320       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL);CHKERRQ(ierr);
321       if (!set_damping_factor_floating) floating_factor = 1.e-12;
322 
323       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",&not_damp_floating,NULL);CHKERRQ(ierr);
324 
325       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",&not_remove_nullspace_floating,NULL);CHKERRQ(ierr);
326 
327       if (pcis->pure_neumann) {  /* floating subdomain */
328         if (!(not_damp_floating)) {
329           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
330           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
331         }
332         if (!(not_remove_nullspace_floating)) {
333           MatNullSpace nullsp;
334           ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr);
335           ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr);
336           ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
337         }
338       } else {  /* fixed subdomain */
339         if (damp_fixed) {
340           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
341           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
342         }
343         if (remove_nullspace_fixed) {
344           MatNullSpace nullsp;
345           ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr);
346           ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr);
347           ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
348         }
349       }
350     }
351     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
352     ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr);
353   }
354 
355   PetscFunctionReturn(0);
356 }
357 
358 /* -------------------------------------------------------------------------- */
359 /*
360    PCISDestroy -
361 */
362 #undef __FUNCT__
363 #define __FUNCT__ "PCISDestroy"
364 PetscErrorCode  PCISDestroy(PC pc)
365 {
366   PC_IS          *pcis = (PC_IS*)(pc->data);
367   PetscErrorCode ierr;
368 
369   PetscFunctionBegin;
370   ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr);
371   ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr);
372   ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr);
373   ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr);
374   ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
375   ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
376   ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
377   ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
378   ierr = VecDestroy(&pcis->D);CHKERRQ(ierr);
379   ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr);
380   ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr);
381   ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr);
382   ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr);
383   ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr);
384   ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr);
385   ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr);
386   ierr = VecDestroy(&pcis->vec4_D);CHKERRQ(ierr);
387   ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr);
388   ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr);
389   ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr);
390   ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr);
391   ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr);
392   ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr);
393   ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr);
394   ierr = PetscFree(pcis->work_N);CHKERRQ(ierr);
395   if (pcis->n_neigh > -1) {
396     ierr = ISLocalToGlobalMappingRestoreInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
397   }
398   ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr);
399   ierr = ISLocalToGlobalMappingDestroy(&pcis->BtoNmap);CHKERRQ(ierr);
400   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL);CHKERRQ(ierr);
401   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL);CHKERRQ(ierr);
402   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL);CHKERRQ(ierr);
403   PetscFunctionReturn(0);
404 }
405 
406 /* -------------------------------------------------------------------------- */
407 /*
408    PCISCreate -
409 */
410 #undef __FUNCT__
411 #define __FUNCT__ "PCISCreate"
412 PetscErrorCode  PCISCreate(PC pc)
413 {
414   PC_IS          *pcis = (PC_IS*)(pc->data);
415   PetscErrorCode ierr;
416 
417   PetscFunctionBegin;
418   pcis->is_B_local       = 0;
419   pcis->is_I_local       = 0;
420   pcis->is_B_global      = 0;
421   pcis->is_I_global      = 0;
422   pcis->A_II             = 0;
423   pcis->A_IB             = 0;
424   pcis->A_BI             = 0;
425   pcis->A_BB             = 0;
426   pcis->D                = 0;
427   pcis->ksp_N            = 0;
428   pcis->ksp_D            = 0;
429   pcis->vec1_N           = 0;
430   pcis->vec2_N           = 0;
431   pcis->vec1_D           = 0;
432   pcis->vec2_D           = 0;
433   pcis->vec3_D           = 0;
434   pcis->vec1_B           = 0;
435   pcis->vec2_B           = 0;
436   pcis->vec3_B           = 0;
437   pcis->vec1_global      = 0;
438   pcis->work_N           = 0;
439   pcis->global_to_D      = 0;
440   pcis->N_to_B           = 0;
441   pcis->global_to_B      = 0;
442   pcis->mapping          = 0;
443   pcis->BtoNmap          = 0;
444   pcis->n_neigh          = -1;
445   pcis->scaling_factor   = 1.0;
446   pcis->reusesubmatrices = PETSC_TRUE;
447   /* composing functions */
448   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS);CHKERRQ(ierr);
449   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr);
450   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr);
451   PetscFunctionReturn(0);
452 }
453 
454 /* -------------------------------------------------------------------------- */
455 /*
456    PCISApplySchur -
457 
458    Input parameters:
459 .  pc - preconditioner context
460 .  v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
461 
462    Output parameters:
463 .  vec1_B - result of Schur complement applied to chunk
464 .  vec2_B - garbage (used as work space), or null (and v is used as workspace)
465 .  vec1_D - garbage (used as work space)
466 .  vec2_D - garbage (used as work space)
467 
468 */
469 #undef __FUNCT__
470 #define __FUNCT__ "PCISApplySchur"
471 PetscErrorCode  PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
472 {
473   PetscErrorCode ierr;
474   PC_IS          *pcis = (PC_IS*)(pc->data);
475 
476   PetscFunctionBegin;
477   if (!vec2_B) vec2_B = v;
478 
479   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
480   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
481   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
482   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
483   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
484   PetscFunctionReturn(0);
485 }
486 
487 /* -------------------------------------------------------------------------- */
488 /*
489    PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
490    including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
491    mode.
492 
493    Input parameters:
494 .  pc - preconditioner context
495 .  array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
496 .  v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
497 
498    Output parameter:
499 .  array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
500 .  v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
501 
502    Notes:
503    The entries in the array that do not correspond to interface nodes remain unaltered.
504 */
505 #undef __FUNCT__
506 #define __FUNCT__ "PCISScatterArrayNToVecB"
507 PetscErrorCode  PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc)
508 {
509   PetscInt       i;
510   const PetscInt *idex;
511   PetscErrorCode ierr;
512   PetscScalar    *array_B;
513   PC_IS          *pcis = (PC_IS*)(pc->data);
514 
515   PetscFunctionBegin;
516   ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr);
517   ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
518 
519   if (smode == SCATTER_FORWARD) {
520     if (imode == INSERT_VALUES) {
521       for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]];
522     } else {  /* ADD_VALUES */
523       for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]];
524     }
525   } else {  /* SCATTER_REVERSE */
526     if (imode == INSERT_VALUES) {
527       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i];
528     } else {  /* ADD_VALUES */
529       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i];
530     }
531   }
532   ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
533   ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr);
534   PetscFunctionReturn(0);
535 }
536 
537 /* -------------------------------------------------------------------------- */
538 /*
539    PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
540    More precisely, solves the problem:
541                                         [ A_II  A_IB ] [ . ]   [ 0 ]
542                                         [            ] [   ] = [   ]
543                                         [ A_BI  A_BB ] [ x ]   [ b ]
544 
545    Input parameters:
546 .  pc - preconditioner context
547 .  b - vector of local interface nodes (including ghosts)
548 
549    Output parameters:
550 .  x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
551        complement to b
552 .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
553 .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
554 
555 */
556 #undef __FUNCT__
557 #define __FUNCT__ "PCISApplyInvSchur"
558 PetscErrorCode  PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
559 {
560   PetscErrorCode ierr;
561   PC_IS          *pcis = (PC_IS*)(pc->data);
562 
563   PetscFunctionBegin;
564   /*
565     Neumann solvers.
566     Applying the inverse of the local Schur complement, i.e, solving a Neumann
567     Problem with zero at the interior nodes of the RHS and extracting the interface
568     part of the solution. inverse Schur complement is applied to b and the result
569     is stored in x.
570   */
571   /* Setting the RHS vec1_N */
572   ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr);
573   ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
574   ierr = VecScatterEnd  (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
575   /* Checking for consistency of the RHS */
576   {
577     PetscBool flg = PETSC_FALSE;
578     ierr = PetscOptionsGetBool(NULL,"-pc_is_check_consistency",&flg,NULL);CHKERRQ(ierr);
579     if (flg) {
580       PetscScalar average;
581       PetscViewer viewer;
582       ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr);
583 
584       ierr    = VecSum(vec1_N,&average);CHKERRQ(ierr);
585       average = average / ((PetscReal)pcis->n);
586       ierr    = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
587       if (pcis->pure_neumann) {
588         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
589       } else {
590         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed.    Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
591       }
592       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
593       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
594     }
595   }
596   /* Solving the system for vec2_N */
597   ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr);
598   /* Extracting the local interface vector out of the solution */
599   ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
600   ierr = VecScatterEnd  (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
601   PetscFunctionReturn(0);
602 }
603 
604 
605 
606 
607 
608 
609 
610 
611 
612