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