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