Actual source code: pcis.c
2: #include src/ksp/pc/impls/is/pcis.h
4: /* -------------------------------------------------------------------------- */
5: /*
6: PCISSetUp -
7: */
10: PetscErrorCode PCISSetUp(PC pc)
11: {
12: PC_IS *pcis = (PC_IS*)(pc->data);
13: Mat_IS *matis = (Mat_IS*)pc->mat->data;
14: PetscInt i;
15: PetscErrorCode ierr;
16: PetscTruth flg;
17:
19: PetscTypeCompare((PetscObject)pc->mat,MATIS,&flg);
20: if (!flg){
21: SETERRQ(PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS");
22: }
24: pcis->pure_neumann = matis->pure_neumann;
26: /*
27: Creating the local vector vec1_N, containing the inverse of the number
28: of subdomains to which each local node (either owned or ghost)
29: pertains. To accomplish that, we scatter local vectors of 1's to
30: a global vector (adding the values); scatter the result back to
31: local vectors and finally invert the result.
32: */
33: {
34: Vec counter;
35: PetscScalar one=1.0, zero=0.0;
36: VecDuplicate(matis->x,&pcis->vec1_N);
37: MatGetVecs(pc->pmat,&counter,0); /* temporary auxiliar vector */
38: VecSet(&zero,counter);
39: VecSet(&one,pcis->vec1_N);
40: VecScatterBegin(pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE,matis->ctx);
41: VecScatterEnd (pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE,matis->ctx);
42: VecScatterBegin(counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD,matis->ctx);
43: VecScatterEnd (counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD,matis->ctx);
44: VecDestroy(counter);
45: }
46: /*
47: Creating local and global index sets for interior and
48: inteface nodes. Notice that interior nodes have D[i]==1.0.
49: */
50: {
51: PetscInt n_I;
52: PetscInt *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global;
53: PetscScalar *array;
54: /* Identifying interior and interface nodes, in local numbering */
55: VecGetSize(pcis->vec1_N,&pcis->n);
56: VecGetArray(pcis->vec1_N,&array);
57: PetscMalloc(pcis->n*sizeof(PetscInt),&idx_I_local);
58: PetscMalloc(pcis->n*sizeof(PetscInt),&idx_B_local);
59: for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) {
60: if (array[i] == 1.0) { idx_I_local[n_I] = i; n_I++; }
61: else { idx_B_local[pcis->n_B] = i; pcis->n_B++; }
62: }
63: /* Getting the global numbering */
64: idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
65: idx_I_global = idx_B_local + pcis->n_B;
66: ISLocalToGlobalMappingApply(matis->mapping,pcis->n_B,idx_B_local,idx_B_global);
67: ISLocalToGlobalMappingApply(matis->mapping,n_I, idx_I_local,idx_I_global);
68: /* Creating the index sets. */
69: ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_local, &pcis->is_B_local);
70: ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_global,&pcis->is_B_global);
71: ISCreateGeneral(MPI_COMM_SELF,n_I ,idx_I_local, &pcis->is_I_local);
72: ISCreateGeneral(MPI_COMM_SELF,n_I ,idx_I_global,&pcis->is_I_global);
73: /* Freeing memory and restoring arrays */
74: PetscFree(idx_B_local);
75: PetscFree(idx_I_local);
76: VecRestoreArray(pcis->vec1_N,&array);
77: }
79: /*
80: Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
81: is such that interior nodes come first than the interface ones, we have
83: [ | ]
84: [ A_II | A_IB ]
85: A = [ | ]
86: [-----------+------]
87: [ A_BI | A_BB ]
88: */
90: MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_II);
91: MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_IB);
92: MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_BI);
93: MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_BB);
95: /*
96: Creating work vectors and arrays
97: */
98: /* pcis->vec1_N has already been created */
99: VecDuplicate(pcis->vec1_N,&pcis->vec2_N);
100: VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);
101: VecDuplicate(pcis->vec1_D,&pcis->vec2_D);
102: VecDuplicate(pcis->vec1_D,&pcis->vec3_D);
103: VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);
104: VecDuplicate(pcis->vec1_B,&pcis->vec2_B);
105: VecDuplicate(pcis->vec1_B,&pcis->vec3_B);
106: MatGetVecs(pc->pmat,&pcis->vec1_global,0);
107: PetscMalloc((pcis->n)*sizeof(PetscScalar),&pcis->work_N);
109: /* Creating the scatter contexts */
110: VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);
111: VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);
112: VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);
114: /* Creating scaling "matrix" D, from information in vec1_N */
115: VecDuplicate(pcis->vec1_B,&pcis->D);
116: VecScatterBegin(pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);
117: VecScatterEnd (pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);
118: VecReciprocal(pcis->D);
120: /* See historical note 01, at the bottom of this file. */
122: /*
123: Creating the KSP contexts for the local Dirichlet and Neumann problems.
124: */
125: {
126: PC pc_ctx;
127: /* Dirichlet */
128: KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);
129: KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II,SAME_PRECONDITIONER);
130: KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");
131: KSPGetPC(pcis->ksp_D,&pc_ctx);
132: PCSetType(pc_ctx,PCLU);
133: KSPSetType(pcis->ksp_D,KSPPREONLY);
134: KSPSetFromOptions(pcis->ksp_D);
135: /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
136: KSPSetUp(pcis->ksp_D);
137: /* Neumann */
138: KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);
139: KSPSetOperators(pcis->ksp_N,matis->A,matis->A,SAME_PRECONDITIONER);
140: KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");
141: KSPGetPC(pcis->ksp_N,&pc_ctx);
142: PCSetType(pc_ctx,PCLU);
143: KSPSetType(pcis->ksp_N,KSPPREONLY);
144: KSPSetFromOptions(pcis->ksp_N);
145: {
146: PetscTruth damp_fixed,
147: remove_nullspace_fixed,
148: set_damping_factor_floating,
149: not_damp_floating,
150: not_remove_nullspace_floating;
151: PetscReal fixed_factor,
152: floating_factor;
154: PetscOptionsGetReal(pc_ctx->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);
155: if (!damp_fixed) { fixed_factor = 0.0; }
156: PetscOptionsHasName(pc_ctx->prefix,"-pc_is_damp_fixed",&damp_fixed);
158: PetscOptionsHasName(pc_ctx->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed);
160: PetscOptionsGetReal(pc_ctx->prefix,"-pc_is_set_damping_factor_floating",
161: &floating_factor,&set_damping_factor_floating);
162: if (!set_damping_factor_floating) { floating_factor = 0.0; }
163: PetscOptionsHasName(pc_ctx->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating);
164: if (!set_damping_factor_floating) { floating_factor = 1.e-12; }
166: PetscOptionsHasName(pc_ctx->prefix,"-pc_is_not_damp_floating",¬_damp_floating);
168: PetscOptionsHasName(pc_ctx->prefix,"-pc_is_not_remove_nullspace_floating",¬_remove_nullspace_floating);
170: if (pcis->pure_neumann) { /* floating subdomain */
171: if (!(not_damp_floating)) {
172: PCLUSetDamping (pc_ctx,floating_factor);
173: PCILUSetDamping(pc_ctx,floating_factor);
174: }
175: if (!(not_remove_nullspace_floating)){
176: MatNullSpace nullsp;
177: MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,PETSC_NULL,&nullsp);
178: KSPSetNullSpace(pcis->ksp_N,nullsp);
179: MatNullSpaceDestroy(nullsp);
180: }
181: } else { /* fixed subdomain */
182: if (damp_fixed) {
183: PCLUSetDamping (pc_ctx,fixed_factor);
184: PCILUSetDamping(pc_ctx,fixed_factor);
185: }
186: if (remove_nullspace_fixed) {
187: MatNullSpace nullsp;
188: MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,PETSC_NULL,&nullsp);
189: KSPSetNullSpace(pcis->ksp_N,nullsp);
190: MatNullSpaceDestroy(nullsp);
191: }
192: }
193: }
194: /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
195: KSPSetUp(pcis->ksp_N);
196: }
198: ISLocalToGlobalMappingGetInfo(((Mat_IS*)(pc->mat->data))->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));
199: pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_TRUE;
200: return(0);
201: }
203: /* -------------------------------------------------------------------------- */
204: /*
205: PCISDestroy -
206: */
209: PetscErrorCode PCISDestroy(PC pc)
210: {
211: PC_IS *pcis = (PC_IS*)(pc->data);
215: if (pcis->is_B_local) {ISDestroy(pcis->is_B_local);}
216: if (pcis->is_I_local) {ISDestroy(pcis->is_I_local);}
217: if (pcis->is_B_global) {ISDestroy(pcis->is_B_global);}
218: if (pcis->is_I_global) {ISDestroy(pcis->is_I_global);}
219: if (pcis->A_II) {MatDestroy(pcis->A_II);}
220: if (pcis->A_IB) {MatDestroy(pcis->A_IB);}
221: if (pcis->A_BI) {MatDestroy(pcis->A_BI);}
222: if (pcis->A_BB) {MatDestroy(pcis->A_BB);}
223: if (pcis->D) {VecDestroy(pcis->D);}
224: if (pcis->ksp_N) {KSPDestroy(pcis->ksp_N);}
225: if (pcis->ksp_D) {KSPDestroy(pcis->ksp_D);}
226: if (pcis->vec1_N) {VecDestroy(pcis->vec1_N);}
227: if (pcis->vec2_N) {VecDestroy(pcis->vec2_N);}
228: if (pcis->vec1_D) {VecDestroy(pcis->vec1_D);}
229: if (pcis->vec2_D) {VecDestroy(pcis->vec2_D);}
230: if (pcis->vec3_D) {VecDestroy(pcis->vec3_D);}
231: if (pcis->vec1_B) {VecDestroy(pcis->vec1_B);}
232: if (pcis->vec2_B) {VecDestroy(pcis->vec2_B);}
233: if (pcis->vec3_B) {VecDestroy(pcis->vec3_B);}
234: if (pcis->vec1_global) {VecDestroy(pcis->vec1_global);}
235: if (pcis->work_N) {PetscFree(pcis->work_N);}
236: if (pcis->global_to_D) {VecScatterDestroy(pcis->global_to_D);}
237: if (pcis->N_to_B) {VecScatterDestroy(pcis->N_to_B);}
238: if (pcis->global_to_B) {VecScatterDestroy(pcis->global_to_B);}
239: if (pcis->ISLocalToGlobalMappingGetInfoWasCalled) {
240: ISLocalToGlobalMappingRestoreInfo((ISLocalToGlobalMapping)0,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));
241: }
242: return(0);
243: }
245: /* -------------------------------------------------------------------------- */
246: /*
247: PCISCreate -
248: */
251: PetscErrorCode PCISCreate(PC pc)
252: {
253: PC_IS *pcis = (PC_IS*)(pc->data);
256: pcis->is_B_local = 0;
257: pcis->is_I_local = 0;
258: pcis->is_B_global = 0;
259: pcis->is_I_global = 0;
260: pcis->A_II = 0;
261: pcis->A_IB = 0;
262: pcis->A_BI = 0;
263: pcis->A_BB = 0;
264: pcis->D = 0;
265: pcis->ksp_N = 0;
266: pcis->ksp_D = 0;
267: pcis->vec1_N = 0;
268: pcis->vec2_N = 0;
269: pcis->vec1_D = 0;
270: pcis->vec2_D = 0;
271: pcis->vec3_D = 0;
272: pcis->vec1_B = 0;
273: pcis->vec2_B = 0;
274: pcis->vec3_B = 0;
275: pcis->vec1_global = 0;
276: pcis->work_N = 0;
277: pcis->global_to_D = 0;
278: pcis->N_to_B = 0;
279: pcis->global_to_B = 0;
280: pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_FALSE;
281: return(0);
282: }
284: /* -------------------------------------------------------------------------- */
285: /*
286: PCISApplySchur -
288: Input parameters:
289: . pc - preconditioner context
290: . v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
292: Output parameters:
293: . vec1_B - result of Schur complement applied to chunk
294: . vec2_B - garbage (used as work space), or null (and v is used as workspace)
295: . vec1_D - garbage (used as work space)
296: . vec2_D - garbage (used as work space)
298: */
301: PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
302: {
304: PetscScalar m_one = -1.0;
305: PC_IS *pcis = (PC_IS*)(pc->data);
308: if (!vec2_B) { vec2_B = v; }
310: MatMult(pcis->A_BB,v,vec1_B);
311: MatMult(pcis->A_IB,v,vec1_D);
312: KSPSolve(pcis->ksp_D,vec1_D,vec2_D);
313: MatMult(pcis->A_BI,vec2_D,vec2_B);
314: VecAXPY(&m_one,vec2_B,vec1_B);
315: return(0);
316: }
318: /* -------------------------------------------------------------------------- */
319: /*
320: PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
321: including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
322: mode.
324: Input parameters:
325: . pc - preconditioner context
326: . array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
327: . v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
329: Output parameter:
330: . array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
331: . v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
333: Notes:
334: The entries in the array that do not correspond to interface nodes remain unaltered.
335: */
338: PetscErrorCode PCISScatterArrayNToVecB (PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc)
339: {
340: PetscInt i, *idex;
342: PetscScalar *array_B;
343: PC_IS *pcis = (PC_IS*)(pc->data);
346: VecGetArray(v_B,&array_B);
347: ISGetIndices(pcis->is_B_local,&idex);
349: if (smode == SCATTER_FORWARD) {
350: if (imode == INSERT_VALUES) {
351: for (i=0; i<pcis->n_B; i++) { array_B[i] = array_N[idex[i]]; }
352: } else { /* ADD_VALUES */
353: for (i=0; i<pcis->n_B; i++) { array_B[i] += array_N[idex[i]]; }
354: }
355: } else { /* SCATTER_REVERSE */
356: if (imode == INSERT_VALUES) {
357: for (i=0; i<pcis->n_B; i++) { array_N[idex[i]] = array_B[i]; }
358: } else { /* ADD_VALUES */
359: for (i=0; i<pcis->n_B; i++) { array_N[idex[i]] += array_B[i]; }
360: }
361: }
362: ISRestoreIndices(pcis->is_B_local,&idex);
363: VecRestoreArray(v_B,&array_B);
364: return(0);
365: }
367: /* -------------------------------------------------------------------------- */
368: /*
369: PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
370: More precisely, solves the problem:
371: [ A_II A_IB ] [ . ] [ 0 ]
372: [ ] [ ] = [ ]
373: [ A_BI A_BB ] [ x ] [ b ]
375: Input parameters:
376: . pc - preconditioner context
377: . b - vector of local interface nodes (including ghosts)
379: Output parameters:
380: . x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
381: complement to b
382: . vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
383: . vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
385: */
388: PetscErrorCode PCISApplyInvSchur (PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
389: {
391: PC_IS *pcis = (PC_IS*)(pc->data);
392: PetscScalar zero = 0.0;
395: /*
396: Neumann solvers.
397: Applying the inverse of the local Schur complement, i.e, solving a Neumann
398: Problem with zero at the interior nodes of the RHS and extracting the interface
399: part of the solution. inverse Schur complement is applied to b and the result
400: is stored in x.
401: */
402: /* Setting the RHS vec1_N */
403: VecSet(&zero,vec1_N);
404: VecScatterBegin(b,vec1_N,INSERT_VALUES,SCATTER_REVERSE,pcis->N_to_B);
405: VecScatterEnd (b,vec1_N,INSERT_VALUES,SCATTER_REVERSE,pcis->N_to_B);
406: /* Checking for consistency of the RHS */
407: {
408: PetscTruth flg;
409: PetscOptionsHasName(PETSC_NULL,"-pc_is_check_consistency",&flg);
410: if (flg) {
411: PetscScalar average;
412: VecSum(vec1_N,&average);
413: average = average / ((PetscReal)pcis->n);
414: if (pcis->pure_neumann) {
415: PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_(pc->comm),"Subdomain %04d is floating. Average = % 1.14e\n",
416: PetscGlobalRank,PetscAbsScalar(average));
417: } else {
418: PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_(pc->comm),"Subdomain %04d is fixed. Average = % 1.14e\n",
419: PetscGlobalRank,PetscAbsScalar(average));
420: }
421: PetscViewerFlush(PETSC_VIEWER_STDOUT_(pc->comm));
422: }
423: }
424: /* Solving the system for vec2_N */
425: KSPSolve(pcis->ksp_N,vec1_N,vec2_N);
426: /* Extracting the local interface vector out of the solution */
427: VecScatterBegin(vec2_N,x,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);
428: VecScatterEnd (vec2_N,x,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);
429: return(0);
430: }