Actual source code: stagmulti.c
1: /* Internal and DMStag-specific functions related to multigrid */
2: #include <petsc/private/dmstagimpl.h>
4: /*@C
5: DMStagRestrictSimple - restricts data from a fine to a coarse DMStag, in the simplest way
7: Values on coarse cells are averages of all fine cells that they cover.
8: Thus, values on vertices are injected, values on edges are averages
9: of the underlying two fine edges, and and values on elements in
10: d dimensions are averages of $2^d$ underlying elements.
12: Input Parameters:
13: + dmf - fine DM
14: . xf - data on fine DM
15: - dmc - coarse DM
17: Output Parameter:
18: . xc - data on coarse DM
20: Level: advanced
22: .seealso: DMRestrict(), DMCoarsen(), DMSTAG, DMCreateInjection()
24: @*/
25: PetscErrorCode DMStagRestrictSimple(DM dmf, Vec xf, DM dmc, Vec xc)
26: {
27: PetscInt dim;
29: DMGetDimension(dmf, &dim);
30: switch (dim) {
31: case 1:
32: DMStagRestrictSimple_1d(dmf, xf, dmc, xc);
33: break;
34: case 2:
35: DMStagRestrictSimple_2d(dmf, xf, dmc, xc);
36: break;
37: case 3:
38: DMStagRestrictSimple_3d(dmf, xf, dmc, xc);
39: break;
40: default:
41: SETERRQ(PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %" PetscInt_FMT "", dim);
42: break;
43: }
44: return 0;
45: }
47: /* Code duplication note: the next two functions are nearly identical, save the inclusion of the element terms */
48: PETSC_INTERN PetscErrorCode DMStagPopulateInterpolation1d_a_b_Private(DM dmc, DM dmf, Mat A)
49: {
50: PetscInt exf, startexf, nexf, nextraxf, startexc;
51: PetscInt dof[2];
52: const PetscInt dim = 1;
54: DMStagGetDOF(dmc, &dof[0], &dof[1], NULL, NULL);
58: /* In 1D, each fine point can receive data from at most 2 coarse points, at most one of which could be off-process */
59: MatSeqAIJSetPreallocation(A, 2, NULL);
60: MatMPIAIJSetPreallocation(A, 2, NULL, 1, NULL);
62: DMStagGetCorners(dmf, &startexf, NULL, NULL, &nexf, NULL, NULL, &nextraxf, NULL, NULL);
63: DMStagGetCorners(dmc, &startexc, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);
64: for (exf = startexf; exf < startexf + nexf + nextraxf; ++exf) {
65: PetscInt exc, exf_local;
66: exf_local = exf - startexf;
67: exc = startexc + exf_local / 2;
68: /* "even" vertices are just injected, odd vertices averaged */
69: if (exf_local % 2 == 0) {
70: DMStagStencil rowf, colc;
71: PetscInt ir, ic;
72: const PetscScalar one = 1.0;
74: rowf.i = exf;
75: rowf.c = 0;
76: rowf.loc = DMSTAG_LEFT;
77: colc.i = exc;
78: colc.c = 0;
79: colc.loc = DMSTAG_LEFT;
80: DMStagStencilToIndexLocal(dmf, dim, 1, &rowf, &ir);
81: DMStagStencilToIndexLocal(dmc, dim, 1, &colc, &ic);
82: MatSetValuesLocal(A, 1, &ir, 1, &ic, &one, INSERT_VALUES);
83: } else {
84: DMStagStencil rowf, colc[2];
85: PetscInt ir, ic[2];
86: const PetscScalar weight[2] = {0.5, 0.5};
88: rowf.i = exf;
89: rowf.c = 0;
90: rowf.loc = DMSTAG_LEFT;
91: colc[0].i = exc;
92: colc[0].c = 0;
93: colc[0].loc = DMSTAG_LEFT;
94: colc[1].i = exc;
95: colc[1].c = 0;
96: colc[1].loc = DMSTAG_RIGHT;
97: DMStagStencilToIndexLocal(dmf, dim, 1, &rowf, &ir);
98: DMStagStencilToIndexLocal(dmc, dim, 2, colc, ic);
99: MatSetValuesLocal(A, 1, &ir, 2, ic, weight, INSERT_VALUES);
100: }
101: /* Elements (excluding "extra" dummies) */
102: if (dof[1] > 0 && exf < startexf + nexf) {
103: DMStagStencil rowf, colc;
104: PetscInt ir, ic;
105: const PetscScalar weight = 1.0;
107: rowf.i = exf;
108: rowf.c = 0;
109: rowf.loc = DMSTAG_ELEMENT; /* Note that this assumes only 1 dof */
110: colc.i = exc;
111: colc.c = 0;
112: colc.loc = DMSTAG_ELEMENT;
113: DMStagStencilToIndexLocal(dmf, dim, 1, &rowf, &ir);
114: DMStagStencilToIndexLocal(dmc, dim, 1, &colc, &ic);
115: MatSetValuesLocal(A, 1, &ir, 1, &ic, &weight, INSERT_VALUES);
116: }
117: }
118: return 0;
119: }
121: PETSC_INTERN PetscErrorCode DMStagPopulateInterpolation2d_0_a_b_Private(DM dmc, DM dmf, Mat A)
122: {
123: PetscInt exf, eyf, startexf, starteyf, nexf, neyf, nextraxf, nextrayf, startexc, starteyc, Nexf, Neyf;
124: PetscInt dof[3];
125: const PetscInt dim = 2;
127: DMStagGetDOF(dmc, &dof[0], &dof[1], &dof[2], NULL);
131: /* In 2D, each fine point can receive data from at most 4 coarse points, at most 3 of which could be off-process */
132: MatSeqAIJSetPreallocation(A, 4, NULL);
133: MatMPIAIJSetPreallocation(A, 4, NULL, 3, NULL);
135: DMStagGetCorners(dmf, &startexf, &starteyf, NULL, &nexf, &neyf, NULL, &nextraxf, &nextrayf, NULL);
136: DMStagGetCorners(dmc, &startexc, &starteyc, NULL, NULL, NULL, NULL, NULL, NULL, NULL);
137: DMStagGetGlobalSizes(dmf, &Nexf, &Neyf, NULL);
138: for (eyf = starteyf; eyf < starteyf + neyf + nextrayf; ++eyf) {
139: PetscInt eyc, eyf_local;
141: eyf_local = eyf - starteyf;
142: eyc = starteyc + eyf_local / 2;
143: for (exf = startexf; exf < startexf + nexf + nextraxf; ++exf) {
144: PetscInt exc, exf_local;
146: exf_local = exf - startexf;
147: exc = startexc + exf_local / 2;
148: /* Left edges (excluding top "extra" dummy row) */
149: if (eyf < starteyf + neyf) {
150: DMStagStencil rowf, colc[4];
151: PetscInt ir, ic[4], nweight;
152: PetscScalar weight[4];
154: rowf.i = exf;
155: rowf.j = eyf;
156: rowf.c = 0;
157: rowf.loc = DMSTAG_LEFT;
158: colc[0].i = exc;
159: colc[0].j = eyc;
160: colc[0].c = 0;
161: colc[0].loc = DMSTAG_LEFT;
162: if (exf_local % 2 == 0) {
163: if (eyf == Neyf - 1 || eyf == 0) {
164: /* Note - this presumes something like a Neumann condition, assuming
165: a ghost edge with the same value as the adjacent physical edge*/
166: nweight = 1;
167: weight[0] = 1.0;
168: } else {
169: nweight = 2;
170: weight[0] = 0.75;
171: weight[1] = 0.25;
172: if (eyf_local % 2 == 0) {
173: colc[1].i = exc;
174: colc[1].j = eyc - 1;
175: colc[1].c = 0;
176: colc[1].loc = DMSTAG_LEFT;
177: } else {
178: colc[1].i = exc;
179: colc[1].j = eyc + 1;
180: colc[1].c = 0;
181: colc[1].loc = DMSTAG_LEFT;
182: }
183: }
184: } else {
185: colc[1].i = exc;
186: colc[1].j = eyc;
187: colc[1].c = 0;
188: colc[1].loc = DMSTAG_RIGHT;
189: if (eyf == Neyf - 1 || eyf == 0) {
190: /* Note - this presumes something like a Neumann condition, assuming
191: a ghost edge with the same value as the adjacent physical edge*/
192: nweight = 2;
193: weight[0] = 0.5;
194: weight[1] = 0.5;
195: } else {
196: nweight = 4;
197: weight[0] = 0.375;
198: weight[1] = 0.375;
199: weight[2] = 0.125;
200: weight[3] = 0.125;
201: if (eyf_local % 2 == 0) {
202: colc[2].i = exc;
203: colc[2].j = eyc - 1;
204: colc[2].c = 0;
205: colc[2].loc = DMSTAG_LEFT;
206: colc[3].i = exc;
207: colc[3].j = eyc - 1;
208: colc[3].c = 0;
209: colc[3].loc = DMSTAG_RIGHT;
210: } else {
211: colc[2].i = exc;
212: colc[2].j = eyc + 1;
213: colc[2].c = 0;
214: colc[2].loc = DMSTAG_LEFT;
215: colc[3].i = exc;
216: colc[3].j = eyc + 1;
217: colc[3].c = 0;
218: colc[3].loc = DMSTAG_RIGHT;
219: }
220: }
221: }
222: DMStagStencilToIndexLocal(dmf, dim, 1, &rowf, &ir);
223: DMStagStencilToIndexLocal(dmc, dim, nweight, colc, ic);
224: MatSetValuesLocal(A, 1, &ir, nweight, ic, weight, INSERT_VALUES);
225: }
226: /* Down edges (excluding right "extra" dummy col) */
227: if (exf < startexf + nexf) {
228: DMStagStencil rowf, colc[4];
229: PetscInt ir, ic[4], nweight;
230: PetscScalar weight[4];
232: rowf.i = exf;
233: rowf.j = eyf;
234: rowf.c = 0;
235: rowf.loc = DMSTAG_DOWN;
236: colc[0].i = exc;
237: colc[0].j = eyc;
238: colc[0].c = 0;
239: colc[0].loc = DMSTAG_DOWN;
240: if (eyf_local % 2 == 0) {
241: if (exf == Nexf - 1 || exf == 0) {
242: /* Note - this presumes something like a Neumann condition, assuming
243: a ghost edge with the same value as the adjacent physical edge*/
244: nweight = 1;
245: weight[0] = 1.0;
246: } else {
247: nweight = 2;
248: weight[0] = 0.75;
249: weight[1] = 0.25;
250: if (exf_local % 2 == 0) {
251: colc[1].i = exc - 1;
252: colc[1].j = eyc;
253: colc[1].c = 0;
254: colc[1].loc = DMSTAG_DOWN;
255: } else {
256: colc[1].i = exc + 1;
257: colc[1].j = eyc;
258: colc[1].c = 0;
259: colc[1].loc = DMSTAG_DOWN;
260: }
261: }
262: } else {
263: colc[1].i = exc;
264: colc[1].j = eyc;
265: colc[1].c = 0;
266: colc[1].loc = DMSTAG_UP;
267: if (exf == Nexf - 1 || exf == 0) {
268: /* Note - this presumes something like a Neumann condition, assuming
269: a ghost edge with the same value as the adjacent physical edge*/
270: nweight = 2;
271: weight[0] = 0.5;
272: weight[1] = 0.5;
273: } else {
274: nweight = 4;
275: weight[0] = 0.375;
276: weight[1] = 0.375;
277: weight[2] = 0.125;
278: weight[3] = 0.125;
279: if (exf_local % 2 == 0) {
280: colc[2].i = exc - 1;
281: colc[2].j = eyc;
282: colc[2].c = 0;
283: colc[2].loc = DMSTAG_DOWN;
284: colc[3].i = exc - 1;
285: colc[3].j = eyc;
286: colc[3].c = 0;
287: colc[3].loc = DMSTAG_UP;
288: } else {
289: colc[2].i = exc + 1;
290: colc[2].j = eyc;
291: colc[2].c = 0;
292: colc[2].loc = DMSTAG_DOWN;
293: colc[3].i = exc + 1;
294: colc[3].j = eyc;
295: colc[3].c = 0;
296: colc[3].loc = DMSTAG_UP;
297: }
298: }
299: }
300: DMStagStencilToIndexLocal(dmf, dim, 1, &rowf, &ir);
301: DMStagStencilToIndexLocal(dmc, dim, nweight, colc, ic);
302: MatSetValuesLocal(A, 1, &ir, nweight, ic, weight, INSERT_VALUES);
303: }
304: /* Elements (excluding "extra" dummy) */
305: if (dof[2] > 0 && exf < startexf + nexf && eyf < starteyf + neyf) {
306: DMStagStencil rowf, colc;
307: PetscInt ir, ic;
308: const PetscScalar weight = 1.0;
310: rowf.i = exf;
311: rowf.j = eyf;
312: rowf.c = 0;
313: rowf.loc = DMSTAG_ELEMENT;
314: colc.i = exc;
315: colc.j = eyc;
316: colc.c = 0;
317: colc.loc = DMSTAG_ELEMENT;
318: DMStagStencilToIndexLocal(dmf, dim, 1, &rowf, &ir);
319: DMStagStencilToIndexLocal(dmc, dim, 1, &colc, &ic);
320: MatSetValuesLocal(A, 1, &ir, 1, &ic, &weight, INSERT_VALUES);
321: }
322: }
323: }
324: return 0;
325: }
327: PETSC_INTERN PetscErrorCode DMStagPopulateInterpolation3d_0_0_a_b_Private(DM dmc, DM dmf, Mat A)
328: {
329: PetscInt exf, eyf, ezf, startexf, starteyf, startezf, nexf, neyf, nezf, nextraxf, nextrayf, nextrazf, startexc, starteyc, startezc, Nexf, Neyf, Nezf;
330: PetscInt dof[4];
331: const PetscInt dim = 3;
334: DMStagGetDOF(dmc, &dof[0], &dof[1], &dof[2], &dof[3]);
338: /* In 3D, each fine point can receive data from at most 8 coarse points, at most 7 of which could be off-process */
339: MatSeqAIJSetPreallocation(A, 8, NULL);
340: MatMPIAIJSetPreallocation(A, 8, NULL, 7, NULL);
342: DMStagGetCorners(dmf, &startexf, &starteyf, &startezf, &nexf, &neyf, &nezf, &nextraxf, &nextrayf, &nextrazf);
343: DMStagGetCorners(dmc, &startexc, &starteyc, &startezc, NULL, NULL, NULL, NULL, NULL, NULL);
344: DMStagGetGlobalSizes(dmf, &Nexf, &Neyf, &Nezf);
345: for (ezf = startezf; ezf < startezf + nezf + nextrayf; ++ezf) {
346: const PetscInt ezf_local = ezf - startezf;
347: const PetscInt ezc = startezc + ezf_local / 2;
348: const PetscBool boundary_z = (PetscBool)(ezf == 0 || ezf == Nezf - 1);
350: for (eyf = starteyf; eyf < starteyf + neyf + nextrayf; ++eyf) {
351: const PetscInt eyf_local = eyf - starteyf;
352: const PetscInt eyc = starteyc + eyf_local / 2;
353: const PetscBool boundary_y = (PetscBool)(eyf == 0 || eyf == Neyf - 1);
355: for (exf = startexf; exf < startexf + nexf + nextraxf; ++exf) {
356: const PetscInt exf_local = exf - startexf;
357: const PetscInt exc = startexc + exf_local / 2;
358: const PetscBool boundary_x = (PetscBool)(exf == 0 || exf == Nexf - 1);
360: /* Left faces (excluding top and front "extra" dummy layers) */
361: if (eyf < starteyf + neyf && ezf < startezf + nezf) {
362: DMStagStencil rowf, colc[8];
363: PetscInt ir, ic[8], nweight;
364: PetscScalar weight[8];
366: rowf.i = exf;
367: rowf.j = eyf;
368: rowf.k = ezf;
369: rowf.c = 0;
370: rowf.loc = DMSTAG_LEFT;
371: colc[0].i = exc;
372: colc[0].j = eyc;
373: colc[0].k = ezc;
374: colc[0].c = 0;
375: colc[0].loc = DMSTAG_LEFT;
376: if (exf_local % 2 == 0) {
377: if (boundary_y) {
378: if (boundary_z) {
379: nweight = 1;
380: weight[0] = 1.0;
381: } else {
382: const PetscInt ezc_offset = ezf_local % 2 == 0 ? -1 : 1;
384: nweight = 2;
385: weight[0] = 0.75;
386: weight[1] = 0.25;
387: colc[1].i = exc;
388: colc[1].j = eyc;
389: colc[1].k = ezc + ezc_offset;
390: colc[1].c = 0;
391: colc[1].loc = DMSTAG_LEFT;
392: }
393: } else {
394: const PetscInt eyc_offset = eyf_local % 2 == 0 ? -1 : 1;
396: if (boundary_z) {
397: nweight = 2;
398: weight[0] = 0.75;
399: weight[1] = 0.25;
400: colc[1].i = exc;
401: colc[1].j = eyc + eyc_offset;
402: colc[1].k = ezc;
403: colc[1].c = 0;
404: colc[1].loc = DMSTAG_LEFT;
405: } else {
406: const PetscInt ezc_offset = ezf_local % 2 == 0 ? -1 : 1;
408: nweight = 4;
409: weight[0] = 0.75 * 0.75;
410: weight[1] = weight[2] = 0.25 * 0.75;
411: weight[3] = 0.25 * 0.25;
412: colc[1].i = exc;
413: colc[1].j = eyc;
414: colc[1].k = ezc + ezc_offset;
415: colc[1].c = 0;
416: colc[1].loc = DMSTAG_LEFT;
417: colc[2].i = exc;
418: colc[2].j = eyc + eyc_offset;
419: colc[2].k = ezc;
420: colc[2].c = 0;
421: colc[2].loc = DMSTAG_LEFT;
422: colc[3].i = exc;
423: colc[3].j = eyc + eyc_offset;
424: colc[3].k = ezc + ezc_offset;
425: colc[3].c = 0;
426: colc[3].loc = DMSTAG_LEFT;
427: }
428: }
429: } else {
430: colc[1].i = exc;
431: colc[1].j = eyc;
432: colc[1].k = ezc;
433: colc[1].c = 0;
434: colc[1].loc = DMSTAG_RIGHT;
435: if (boundary_y) {
436: if (boundary_z) {
437: nweight = 2;
438: weight[0] = weight[1] = 0.5;
439: } else {
440: const PetscInt ezc_offset = ezf_local % 2 == 0 ? -1 : 1;
442: nweight = 4;
443: weight[0] = weight[1] = 0.5 * 0.75;
444: weight[2] = weight[3] = 0.5 * 0.25;
445: colc[2].i = exc;
446: colc[2].j = eyc;
447: colc[2].k = ezc + ezc_offset;
448: colc[2].c = 0;
449: colc[2].loc = DMSTAG_LEFT;
450: colc[3].i = exc;
451: colc[3].j = eyc;
452: colc[3].k = ezc + ezc_offset;
453: colc[3].c = 0;
454: colc[3].loc = DMSTAG_RIGHT;
455: }
456: } else {
457: const PetscInt eyc_offset = eyf_local % 2 == 0 ? -1 : 1;
459: if (boundary_z) {
460: nweight = 4;
461: weight[0] = weight[1] = 0.5 * 0.75;
462: weight[2] = weight[3] = 0.5 * 0.25;
463: colc[2].i = exc;
464: colc[2].j = eyc + eyc_offset;
465: colc[2].k = ezc;
466: colc[2].c = 0;
467: colc[2].loc = DMSTAG_LEFT;
468: colc[3].i = exc;
469: colc[3].j = eyc + eyc_offset;
470: colc[3].k = ezc;
471: colc[3].c = 0;
472: colc[3].loc = DMSTAG_RIGHT;
473: } else {
474: const PetscInt ezc_offset = ezf_local % 2 == 0 ? -1 : 1;
476: nweight = 8;
477: weight[0] = weight[1] = 0.5 * 0.75 * 0.75;
478: weight[2] = weight[3] = weight[4] = weight[5] = 0.5 * 0.25 * 0.75;
479: weight[6] = weight[7] = 0.5 * 0.25 * 0.25;
480: colc[2].i = exc;
481: colc[2].j = eyc;
482: colc[2].k = ezc + ezc_offset;
483: colc[2].c = 0;
484: colc[2].loc = DMSTAG_LEFT;
485: colc[3].i = exc;
486: colc[3].j = eyc;
487: colc[3].k = ezc + ezc_offset;
488: colc[3].c = 0;
489: colc[3].loc = DMSTAG_RIGHT;
490: colc[4].i = exc;
491: colc[4].j = eyc + eyc_offset;
492: colc[4].k = ezc;
493: colc[4].c = 0;
494: colc[4].loc = DMSTAG_LEFT;
495: colc[5].i = exc;
496: colc[5].j = eyc + eyc_offset;
497: colc[5].k = ezc;
498: colc[5].c = 0;
499: colc[5].loc = DMSTAG_RIGHT;
500: colc[6].i = exc;
501: colc[6].j = eyc + eyc_offset;
502: colc[6].k = ezc + ezc_offset;
503: colc[6].c = 0;
504: colc[6].loc = DMSTAG_LEFT;
505: colc[7].i = exc;
506: colc[7].j = eyc + eyc_offset;
507: colc[7].k = ezc + ezc_offset;
508: colc[7].c = 0;
509: colc[7].loc = DMSTAG_RIGHT;
510: }
511: }
512: }
513: DMStagStencilToIndexLocal(dmf, dim, 1, &rowf, &ir);
514: DMStagStencilToIndexLocal(dmc, dim, nweight, colc, ic);
515: MatSetValuesLocal(A, 1, &ir, nweight, ic, weight, INSERT_VALUES);
516: }
518: /* Bottom faces (excluding left and front "extra" dummy layers) */
519: if (exf < startexf + nexf && ezf < startezf + nezf) {
520: DMStagStencil rowf, colc[8];
521: PetscInt ir, ic[8], nweight;
522: PetscScalar weight[8];
524: rowf.i = exf;
525: rowf.j = eyf;
526: rowf.k = ezf;
527: rowf.c = 0;
528: rowf.loc = DMSTAG_DOWN;
529: colc[0].i = exc;
530: colc[0].j = eyc;
531: colc[0].k = ezc;
532: colc[0].c = 0;
533: colc[0].loc = DMSTAG_DOWN;
534: if (eyf_local % 2 == 0) {
535: if (boundary_x) {
536: if (boundary_z) {
537: nweight = 1;
538: weight[0] = 1.0;
539: } else {
540: const PetscInt ezc_offset = ezf_local % 2 == 0 ? -1 : 1;
542: nweight = 2;
543: weight[0] = 0.75;
544: weight[1] = 0.25;
545: colc[1].i = exc;
546: colc[1].j = eyc;
547: colc[1].k = ezc + ezc_offset;
548: colc[1].c = 0;
549: colc[1].loc = DMSTAG_DOWN;
550: }
551: } else {
552: const PetscInt exc_offset = exf_local % 2 == 0 ? -1 : 1;
554: if (boundary_z) {
555: nweight = 2;
556: weight[0] = 0.75;
557: weight[1] = 0.25;
558: colc[1].i = exc + exc_offset;
559: colc[1].j = eyc;
560: colc[1].k = ezc;
561: colc[1].c = 0;
562: colc[1].loc = DMSTAG_DOWN;
563: } else {
564: const PetscInt ezc_offset = ezf_local % 2 == 0 ? -1 : 1;
566: nweight = 4;
567: weight[0] = 0.75 * 0.75;
568: weight[1] = weight[2] = 0.25 * 0.75;
569: weight[3] = 0.25 * 0.25;
570: colc[1].i = exc;
571: colc[1].j = eyc;
572: colc[1].k = ezc + ezc_offset;
573: colc[1].c = 0;
574: colc[1].loc = DMSTAG_DOWN;
575: colc[2].i = exc + exc_offset;
576: colc[2].j = eyc;
577: colc[2].k = ezc;
578: colc[2].c = 0;
579: colc[2].loc = DMSTAG_DOWN;
580: colc[3].i = exc + exc_offset;
581: colc[3].j = eyc;
582: colc[3].k = ezc + ezc_offset;
583: colc[3].c = 0;
584: colc[3].loc = DMSTAG_DOWN;
585: }
586: }
587: } else {
588: colc[1].i = exc;
589: colc[1].j = eyc;
590: colc[1].k = ezc;
591: colc[1].c = 0;
592: colc[1].loc = DMSTAG_UP;
593: if (boundary_x) {
594: if (boundary_z) {
595: nweight = 2;
596: weight[0] = weight[1] = 0.5;
597: } else {
598: const PetscInt ezc_offset = ezf_local % 2 == 0 ? -1 : 1;
600: nweight = 4;
601: weight[0] = weight[1] = 0.5 * 0.75;
602: weight[2] = weight[3] = 0.5 * 0.25;
603: colc[2].i = exc;
604: colc[2].j = eyc;
605: colc[2].k = ezc + ezc_offset;
606: colc[2].c = 0;
607: colc[2].loc = DMSTAG_DOWN;
608: colc[3].i = exc;
609: colc[3].j = eyc;
610: colc[3].k = ezc + ezc_offset;
611: colc[3].c = 0;
612: colc[3].loc = DMSTAG_UP;
613: }
614: } else {
615: const PetscInt exc_offset = exf_local % 2 == 0 ? -1 : 1;
617: if (boundary_z) {
618: nweight = 4;
619: weight[0] = weight[1] = 0.5 * 0.75;
620: weight[2] = weight[3] = 0.5 * 0.25;
621: colc[2].i = exc + exc_offset;
622: colc[2].j = eyc;
623: colc[2].k = ezc;
624: colc[2].c = 0;
625: colc[2].loc = DMSTAG_DOWN;
626: colc[3].i = exc + exc_offset;
627: colc[3].j = eyc;
628: colc[3].k = ezc;
629: colc[3].c = 0;
630: colc[3].loc = DMSTAG_UP;
631: } else {
632: const PetscInt ezc_offset = ezf_local % 2 == 0 ? -1 : 1;
634: nweight = 8;
635: weight[0] = weight[1] = 0.5 * 0.75 * 0.75;
636: weight[2] = weight[3] = weight[4] = weight[5] = 0.5 * 0.25 * 0.75;
637: weight[6] = weight[7] = 0.5 * 0.25 * 0.25;
638: colc[2].i = exc;
639: colc[2].j = eyc;
640: colc[2].k = ezc + ezc_offset;
641: colc[2].c = 0;
642: colc[2].loc = DMSTAG_DOWN;
643: colc[3].i = exc;
644: colc[3].j = eyc;
645: colc[3].k = ezc + ezc_offset;
646: colc[3].c = 0;
647: colc[3].loc = DMSTAG_UP;
648: colc[4].i = exc + exc_offset;
649: colc[4].j = eyc;
650: colc[4].k = ezc;
651: colc[4].c = 0;
652: colc[4].loc = DMSTAG_DOWN;
653: colc[5].i = exc + exc_offset;
654: colc[5].j = eyc;
655: colc[5].k = ezc;
656: colc[5].c = 0;
657: colc[5].loc = DMSTAG_UP;
658: colc[6].i = exc + exc_offset;
659: colc[6].j = eyc;
660: colc[6].k = ezc + ezc_offset;
661: colc[6].c = 0;
662: colc[6].loc = DMSTAG_DOWN;
663: colc[7].i = exc + exc_offset;
664: colc[7].j = eyc;
665: colc[7].k = ezc + ezc_offset;
666: colc[7].c = 0;
667: colc[7].loc = DMSTAG_UP;
668: }
669: }
670: }
671: DMStagStencilToIndexLocal(dmf, dim, 1, &rowf, &ir);
672: DMStagStencilToIndexLocal(dmc, dim, nweight, colc, ic);
673: MatSetValuesLocal(A, 1, &ir, nweight, ic, weight, INSERT_VALUES);
674: }
676: /* Back faces (excluding left and top "extra" dummy layers) */
677: if (exf < startexf + nexf && ezf < startezf + nezf) {
678: DMStagStencil rowf, colc[8];
679: PetscInt ir, ic[8], nweight;
680: PetscScalar weight[8];
682: rowf.i = exf;
683: rowf.j = eyf;
684: rowf.k = ezf;
685: rowf.c = 0;
686: rowf.loc = DMSTAG_BACK;
687: colc[0].i = exc;
688: colc[0].j = eyc;
689: colc[0].k = ezc;
690: colc[0].c = 0;
691: colc[0].loc = DMSTAG_BACK;
692: if (ezf_local % 2 == 0) {
693: if (boundary_x) {
694: if (boundary_y) {
695: nweight = 1;
696: weight[0] = 1.0;
697: } else {
698: const PetscInt eyc_offset = eyf_local % 2 == 0 ? -1 : 1;
700: nweight = 2;
701: weight[0] = 0.75;
702: weight[1] = 0.25;
703: colc[1].i = exc;
704: colc[1].j = eyc + eyc_offset;
705: colc[1].k = ezc;
706: colc[1].c = 0;
707: colc[1].loc = DMSTAG_BACK;
708: }
709: } else {
710: const PetscInt exc_offset = exf_local % 2 == 0 ? -1 : 1;
712: if (boundary_y) {
713: nweight = 2;
714: weight[0] = 0.75;
715: weight[1] = 0.25;
716: colc[1].i = exc + exc_offset;
717: colc[1].j = eyc;
718: colc[1].k = ezc;
719: colc[1].c = 0;
720: colc[1].loc = DMSTAG_BACK;
721: } else {
722: const PetscInt eyc_offset = eyf_local % 2 == 0 ? -1 : 1;
724: nweight = 4;
725: weight[0] = 0.75 * 0.75;
726: weight[1] = weight[2] = 0.25 * 0.75;
727: weight[3] = 0.25 * 0.25;
728: colc[1].i = exc + exc_offset;
729: colc[1].j = eyc;
730: colc[1].k = ezc;
731: colc[1].c = 0;
732: colc[1].loc = DMSTAG_BACK;
733: colc[2].i = exc;
734: colc[2].j = eyc + eyc_offset;
735: colc[2].k = ezc;
736: colc[2].c = 0;
737: colc[2].loc = DMSTAG_BACK;
738: colc[3].i = exc + exc_offset;
739: colc[3].j = eyc + eyc_offset;
740: colc[3].k = ezc;
741: colc[3].c = 0;
742: colc[3].loc = DMSTAG_BACK;
743: }
744: }
745: } else {
746: colc[1].i = exc;
747: colc[1].j = eyc;
748: colc[1].k = ezc;
749: colc[1].c = 0;
750: colc[1].loc = DMSTAG_FRONT;
751: if (boundary_x) {
752: if (boundary_y) {
753: nweight = 2;
754: weight[0] = weight[1] = 0.5;
755: } else {
756: const PetscInt eyc_offset = eyf_local % 2 == 0 ? -1 : 1;
758: nweight = 4;
759: weight[0] = weight[1] = 0.5 * 0.75;
760: weight[2] = weight[3] = 0.5 * 0.25;
761: colc[2].i = exc;
762: colc[2].j = eyc + eyc_offset;
763: colc[2].k = ezc;
764: colc[2].c = 0;
765: colc[2].loc = DMSTAG_BACK;
766: colc[3].i = exc;
767: colc[3].j = eyc + eyc_offset;
768: colc[3].k = ezc;
769: colc[3].c = 0;
770: colc[3].loc = DMSTAG_FRONT;
771: }
772: } else {
773: const PetscInt exc_offset = exf_local % 2 == 0 ? -1 : 1;
775: if (boundary_y) {
776: nweight = 4;
777: weight[0] = weight[1] = 0.5 * 0.75;
778: weight[2] = weight[3] = 0.5 * 0.25;
779: colc[2].i = exc + exc_offset;
780: colc[2].j = eyc;
781: colc[2].k = ezc;
782: colc[2].c = 0;
783: colc[2].loc = DMSTAG_BACK;
784: colc[3].i = exc + exc_offset;
785: colc[3].j = eyc;
786: colc[3].k = ezc;
787: colc[3].c = 0;
788: colc[3].loc = DMSTAG_FRONT;
789: } else {
790: const PetscInt eyc_offset = eyf_local % 2 == 0 ? -1 : 1;
792: nweight = 8;
793: weight[0] = weight[1] = 0.5 * 0.75 * 0.75;
794: weight[2] = weight[3] = weight[4] = weight[5] = 0.5 * 0.25 * 0.75;
795: weight[6] = weight[7] = 0.5 * 0.25 * 0.25;
796: colc[2].i = exc;
797: colc[2].j = eyc + eyc_offset;
798: colc[2].k = ezc;
799: colc[2].c = 0;
800: colc[2].loc = DMSTAG_BACK;
801: colc[3].i = exc;
802: colc[3].j = eyc + eyc_offset;
803: colc[3].k = ezc;
804: colc[3].c = 0;
805: colc[3].loc = DMSTAG_FRONT;
806: colc[4].i = exc + exc_offset;
807: colc[4].j = eyc;
808: colc[4].k = ezc;
809: colc[4].c = 0;
810: colc[4].loc = DMSTAG_BACK;
811: colc[5].i = exc + exc_offset;
812: colc[5].j = eyc;
813: colc[5].k = ezc;
814: colc[5].c = 0;
815: colc[5].loc = DMSTAG_FRONT;
816: colc[6].i = exc + exc_offset;
817: colc[6].j = eyc + eyc_offset;
818: colc[6].k = ezc;
819: colc[6].c = 0;
820: colc[6].loc = DMSTAG_BACK;
821: colc[7].i = exc + exc_offset;
822: colc[7].j = eyc + eyc_offset;
823: colc[7].k = ezc;
824: colc[7].c = 0;
825: colc[7].loc = DMSTAG_FRONT;
826: }
827: }
828: }
829: DMStagStencilToIndexLocal(dmf, dim, 1, &rowf, &ir);
830: DMStagStencilToIndexLocal(dmc, dim, nweight, colc, ic);
831: MatSetValuesLocal(A, 1, &ir, nweight, ic, weight, INSERT_VALUES);
832: }
833: /* Elements */
834: if (dof[3] == 1 && exf < startexf + nexf && eyf < starteyf + neyf && ezf < startezf + nezf) {
835: DMStagStencil rowf, colc;
836: PetscInt ir, ic;
837: const PetscScalar weight = 1.0;
839: rowf.i = exf;
840: rowf.j = eyf;
841: rowf.k = ezf;
842: rowf.c = 0;
843: rowf.loc = DMSTAG_ELEMENT;
844: colc.i = exc;
845: colc.j = eyc;
846: colc.k = ezc;
847: colc.c = 0;
848: colc.loc = DMSTAG_ELEMENT;
849: DMStagStencilToIndexLocal(dmf, dim, 1, &rowf, &ir);
850: DMStagStencilToIndexLocal(dmc, dim, 1, &colc, &ic);
851: MatSetValuesLocal(A, 1, &ir, 1, &ic, &weight, INSERT_VALUES);
852: }
853: }
854: }
855: }
856: return 0;
857: }
859: PETSC_INTERN PetscErrorCode DMStagPopulateRestriction1d_a_b_Private(DM dmc, DM dmf, Mat A)
860: {
861: PetscInt exf, startexf, nexf, nextraxf, startexc, Nexf;
862: PetscInt dof[2];
863: const PetscInt dim = 1;
865: DMStagGetDOF(dmc, &dof[0], &dof[1], NULL, NULL);
869: /* In 1D, each coarse point can receive from up to 3 fine points, one of which may be off-rank */
870: MatSeqAIJSetPreallocation(A, 3, NULL);
871: MatMPIAIJSetPreallocation(A, 3, NULL, 1, NULL);
873: DMStagGetCorners(dmf, &startexf, NULL, NULL, &nexf, NULL, NULL, &nextraxf, NULL, NULL);
874: DMStagGetCorners(dmc, &startexc, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);
875: DMStagGetGlobalSizes(dmf, &Nexf, NULL, NULL);
876: for (exf = startexf; exf < startexf + nexf + nextraxf; ++exf) {
877: PetscInt exc, exf_local;
879: exf_local = exf - startexf;
880: exc = startexc + exf_local / 2;
881: /* "even" vertices contribute to the overlying coarse vertex, odd vertices to the two adjacent */
882: if (exf_local % 2 == 0) {
883: DMStagStencil colf, rowc;
884: PetscInt ir, ic;
885: PetscScalar weight;
887: colf.i = exf;
888: colf.c = 0;
889: colf.loc = DMSTAG_LEFT;
890: rowc.i = exc;
891: rowc.c = 0;
892: rowc.loc = DMSTAG_LEFT;
893: DMStagStencilToIndexLocal(dmc, dim, 1, &rowc, &ir);
894: DMStagStencilToIndexLocal(dmf, dim, 1, &colf, &ic);
895: weight = (exf == Nexf || exf == 0) ? 0.75 : 0.5; /* Assume a Neuman-type condition */
896: MatSetValuesLocal(A, 1, &ir, 1, &ic, &weight, INSERT_VALUES);
897: } else {
898: DMStagStencil colf, rowc[2];
899: PetscInt ic, ir[2];
900: const PetscScalar weight[2] = {0.25, 0.25};
902: colf.i = exf;
903: colf.c = 0;
904: colf.loc = DMSTAG_LEFT;
905: rowc[0].i = exc;
906: rowc[0].c = 0;
907: rowc[0].loc = DMSTAG_LEFT;
908: rowc[1].i = exc;
909: rowc[1].c = 0;
910: rowc[1].loc = DMSTAG_RIGHT;
911: DMStagStencilToIndexLocal(dmc, dim, 2, rowc, ir);
912: DMStagStencilToIndexLocal(dmf, dim, 1, &colf, &ic);
913: MatSetValuesLocal(A, 2, ir, 1, &ic, weight, INSERT_VALUES);
914: }
915: if (dof[1] > 0 && exf < startexf + nexf) {
916: DMStagStencil rowc, colf;
917: PetscInt ir, ic;
918: const PetscScalar weight = 0.5;
920: rowc.i = exc;
921: rowc.c = 0;
922: rowc.loc = DMSTAG_ELEMENT;
923: colf.i = exf;
924: colf.c = 0;
925: colf.loc = DMSTAG_ELEMENT;
926: DMStagStencilToIndexLocal(dmc, dim, 1, &rowc, &ir);
927: DMStagStencilToIndexLocal(dmf, dim, 1, &colf, &ic);
928: MatSetValuesLocal(A, 1, &ir, 1, &ic, &weight, INSERT_VALUES);
929: }
930: }
931: return 0;
932: }
934: PETSC_INTERN PetscErrorCode DMStagPopulateRestriction2d_0_a_b_Private(DM dmc, DM dmf, Mat A)
935: {
936: PetscInt exf, eyf, startexf, starteyf, nexf, neyf, nextraxf, nextrayf, startexc, starteyc, Nexf, Neyf;
937: PetscInt dof[3];
938: const PetscInt dim = 2;
940: DMStagGetDOF(dmc, &dof[0], &dof[1], &dof[2], NULL);
944: /* In 2D, each coarse point can receive from up to 6 fine points,
945: up to 2 of which may be off rank */
946: MatSeqAIJSetPreallocation(A, 6, NULL);
947: MatMPIAIJSetPreallocation(A, 6, NULL, 2, NULL);
949: DMStagGetCorners(dmf, &startexf, &starteyf, NULL, &nexf, &neyf, NULL, &nextraxf, &nextrayf, NULL);
950: DMStagGetCorners(dmc, &startexc, &starteyc, NULL, NULL, NULL, NULL, NULL, NULL, NULL);
951: DMStagGetGlobalSizes(dmf, &Nexf, &Neyf, NULL);
952: for (eyf = starteyf; eyf < starteyf + neyf + nextrayf; ++eyf) {
953: PetscInt eyc, eyf_local;
954: eyf_local = eyf - starteyf;
955: eyc = starteyc + eyf_local / 2;
956: for (exf = startexf; exf < startexf + nexf + nextraxf; ++exf) {
957: PetscInt exc, exf_local;
958: exf_local = exf - startexf;
959: exc = startexc + exf_local / 2;
960: /* Left edges (excluding top "extra" dummy row) */
961: if (eyf < starteyf + neyf) {
962: DMStagStencil rowc[2], colf;
963: PetscInt ir[2], ic, nweight;
964: PetscScalar weight[2];
966: colf.i = exf;
967: colf.j = eyf;
968: colf.c = 0;
969: colf.loc = DMSTAG_LEFT;
970: rowc[0].i = exc;
971: rowc[0].j = eyc;
972: rowc[0].c = 0;
973: rowc[0].loc = DMSTAG_LEFT;
974: if (exf_local % 2 == 0) {
975: nweight = 1;
976: if (exf == Nexf || exf == 0) {
977: /* Note - this presumes something like a Neumann condition, assuming
978: a ghost edge with the same value as the adjacent physical edge*/
979: weight[0] = 0.375;
980: } else {
981: weight[0] = 0.25;
982: }
983: } else {
984: nweight = 2;
985: rowc[1].i = exc;
986: rowc[1].j = eyc;
987: rowc[1].c = 0;
988: rowc[1].loc = DMSTAG_RIGHT;
989: weight[0] = 0.125;
990: weight[1] = 0.125;
991: }
992: DMStagStencilToIndexLocal(dmc, dim, nweight, rowc, ir);
993: DMStagStencilToIndexLocal(dmf, dim, 1, &colf, &ic);
994: MatSetValuesLocal(A, nweight, ir, 1, &ic, weight, INSERT_VALUES);
995: }
996: /* Down edges (excluding right "extra" dummy col) */
997: if (exf < startexf + nexf) {
998: DMStagStencil rowc[2], colf;
999: PetscInt ir[2], ic, nweight;
1000: PetscScalar weight[2];
1002: colf.i = exf;
1003: colf.j = eyf;
1004: colf.c = 0;
1005: colf.loc = DMSTAG_DOWN;
1006: rowc[0].i = exc;
1007: rowc[0].j = eyc;
1008: rowc[0].c = 0;
1009: rowc[0].loc = DMSTAG_DOWN;
1010: if (eyf_local % 2 == 0) {
1011: nweight = 1;
1012: if (eyf == Neyf || eyf == 0) {
1013: /* Note - this presumes something like a Neumann condition, assuming
1014: a ghost edge with the same value as the adjacent physical edge*/
1015: weight[0] = 0.375;
1016: } else {
1017: weight[0] = 0.25;
1018: }
1019: } else {
1020: nweight = 2;
1021: rowc[1].i = exc;
1022: rowc[1].j = eyc;
1023: rowc[1].c = 0;
1024: rowc[1].loc = DMSTAG_UP;
1025: weight[0] = 0.125;
1026: weight[1] = 0.125;
1027: }
1028: DMStagStencilToIndexLocal(dmc, dim, nweight, rowc, ir);
1029: DMStagStencilToIndexLocal(dmf, dim, 1, &colf, &ic);
1030: MatSetValuesLocal(A, nweight, ir, 1, &ic, weight, INSERT_VALUES);
1031: }
1032: /* Elements (excluding "extra" dummies) */
1033: if (dof[2] > 0 && exf < startexf + nexf && eyf < starteyf + neyf) {
1034: DMStagStencil rowc, colf;
1035: PetscInt ir, ic;
1036: const PetscScalar cellScale = 0.25;
1038: rowc.i = exc;
1039: rowc.j = eyc;
1040: rowc.c = 0;
1041: rowc.loc = DMSTAG_ELEMENT;
1042: colf.i = exf;
1043: colf.j = eyf;
1044: colf.c = 0;
1045: colf.loc = DMSTAG_ELEMENT;
1046: DMStagStencilToIndexLocal(dmc, dim, 1, &rowc, &ir);
1047: DMStagStencilToIndexLocal(dmf, dim, 1, &colf, &ic);
1048: MatSetValuesLocal(A, 1, &ir, 1, &ic, &cellScale, INSERT_VALUES);
1049: }
1050: }
1051: }
1052: return 0;
1053: }
1055: PETSC_INTERN PetscErrorCode DMStagPopulateRestriction3d_0_0_a_b_Private(DM dmc, DM dmf, Mat A)
1056: {
1057: PetscInt exf, eyf, ezf, startexf, starteyf, startezf, nexf, neyf, nezf, nextraxf, nextrayf, nextrazf, startexc, starteyc, startezc, Nexf, Neyf, Nezf;
1058: PetscInt dof[4];
1059: const PetscInt dim = 3;
1062: DMStagGetDOF(dmc, &dof[0], &dof[1], &dof[2], &dof[3]);
1066: /* In 3D, each coarse point can receive from up to 12 fine points,
1067: up to 4 of which may be off rank */
1068: MatSeqAIJSetPreallocation(A, 12, NULL);
1069: MatMPIAIJSetPreallocation(A, 12, NULL, 4, NULL);
1071: DMStagGetCorners(dmf, &startexf, &starteyf, &startezf, &nexf, &neyf, &nezf, &nextraxf, &nextrayf, &nextrazf);
1072: DMStagGetCorners(dmc, &startexc, &starteyc, &startezc, NULL, NULL, NULL, NULL, NULL, NULL);
1073: DMStagGetGlobalSizes(dmf, &Nexf, &Neyf, &Nezf);
1075: for (ezf = startezf; ezf < startezf + nezf + nextrazf; ++ezf) {
1076: const PetscInt ezf_local = ezf - startezf;
1077: const PetscInt ezc = startezc + ezf_local / 2;
1079: for (eyf = starteyf; eyf < starteyf + neyf + nextrayf; ++eyf) {
1080: const PetscInt eyf_local = eyf - starteyf;
1081: const PetscInt eyc = starteyc + eyf_local / 2;
1083: for (exf = startexf; exf < startexf + nexf + nextraxf; ++exf) {
1084: const PetscInt exf_local = exf - startexf;
1085: const PetscInt exc = startexc + exf_local / 2;
1087: /* Left faces (excluding top and front "extra" dummy layers) */
1088: if (eyf < starteyf + neyf && ezf < startezf + nezf) {
1089: DMStagStencil rowc[2], colf;
1090: PetscInt ir[2], ic, nweight;
1091: PetscScalar weight[2];
1093: colf.i = exf;
1094: colf.j = eyf;
1095: colf.k = ezf;
1096: colf.c = 0;
1097: colf.loc = DMSTAG_LEFT;
1098: rowc[0].i = exc;
1099: rowc[0].j = eyc;
1100: rowc[0].k = ezc;
1101: rowc[0].c = 0;
1102: rowc[0].loc = DMSTAG_LEFT;
1103: if (exf_local % 2 == 0) {
1104: nweight = 1;
1105: if (exf == Nexf || exf == 0) {
1106: /* Note - this presumes something like a Neumann condition, assuming
1107: a ghost edge with the same value as the adjacent physical edge*/
1108: weight[0] = 0.1875;
1109: } else {
1110: weight[0] = 0.125;
1111: }
1112: } else {
1113: nweight = 2;
1114: rowc[1].i = exc;
1115: rowc[1].j = eyc;
1116: rowc[1].k = ezc;
1117: rowc[1].c = 0;
1118: rowc[1].loc = DMSTAG_RIGHT;
1119: weight[0] = 0.0625;
1120: weight[1] = 0.0625;
1121: }
1122: DMStagStencilToIndexLocal(dmc, dim, nweight, rowc, ir);
1123: DMStagStencilToIndexLocal(dmf, dim, 1, &colf, &ic);
1124: MatSetValuesLocal(A, nweight, ir, 1, &ic, weight, INSERT_VALUES);
1125: }
1127: /* Down faces (excluding right and front "extra" dummy layers) */
1128: if (exf < startexf + nexf && ezf < startezf + nezf) {
1129: DMStagStencil rowc[2], colf;
1130: PetscInt ir[2], ic, nweight;
1131: PetscScalar weight[2];
1133: colf.i = exf;
1134: colf.j = eyf;
1135: colf.k = ezf;
1136: colf.c = 0;
1137: colf.loc = DMSTAG_DOWN;
1138: rowc[0].i = exc;
1139: rowc[0].j = eyc;
1140: rowc[0].k = ezc;
1141: rowc[0].c = 0;
1142: rowc[0].loc = DMSTAG_DOWN;
1143: if (eyf_local % 2 == 0) {
1144: nweight = 1;
1145: if (eyf == Neyf || eyf == 0) {
1146: /* Note - this presumes something like a Neumann condition, assuming
1147: a ghost edge with the same value as the adjacent physical edge*/
1148: weight[0] = 0.1875;
1149: } else {
1150: weight[0] = 0.125;
1151: }
1152: } else {
1153: nweight = 2;
1154: rowc[1].i = exc;
1155: rowc[1].j = eyc;
1156: rowc[1].k = ezc;
1157: rowc[1].c = 0;
1158: rowc[1].loc = DMSTAG_UP;
1159: weight[0] = 0.0625;
1160: weight[1] = 0.0625;
1161: }
1162: DMStagStencilToIndexLocal(dmc, dim, nweight, rowc, ir);
1163: DMStagStencilToIndexLocal(dmf, dim, 1, &colf, &ic);
1164: MatSetValuesLocal(A, nweight, ir, 1, &ic, weight, INSERT_VALUES);
1165: }
1167: /* Back faces (excluding left and top "extra" dummy laers) */
1168: if (exf < startexf + nexf && eyf < starteyf + neyf) {
1169: DMStagStencil rowc[2], colf;
1170: PetscInt ir[2], ic, nweight;
1171: PetscScalar weight[2];
1173: colf.i = exf;
1174: colf.j = eyf;
1175: colf.k = ezf;
1176: colf.c = 0;
1177: colf.loc = DMSTAG_BACK;
1178: rowc[0].i = exc;
1179: rowc[0].j = eyc;
1180: rowc[0].k = ezc;
1181: rowc[0].c = 0;
1182: rowc[0].loc = DMSTAG_BACK;
1183: if (ezf_local % 2 == 0) {
1184: nweight = 1;
1185: if (ezf == Nezf || ezf == 0) {
1186: /* Note - this presumes something like a Neumann condition, assuming
1187: a ghost edge with the same value as the adjacent physical edge*/
1188: weight[0] = 0.1875;
1189: } else {
1190: weight[0] = 0.125;
1191: }
1192: } else {
1193: nweight = 2;
1194: rowc[1].i = exc;
1195: rowc[1].j = eyc;
1196: rowc[1].k = ezc;
1197: rowc[1].c = 0;
1198: rowc[1].loc = DMSTAG_FRONT;
1199: weight[0] = 0.0625;
1200: weight[1] = 0.0625;
1201: }
1202: DMStagStencilToIndexLocal(dmc, dim, nweight, rowc, ir);
1203: DMStagStencilToIndexLocal(dmf, dim, 1, &colf, &ic);
1204: MatSetValuesLocal(A, nweight, ir, 1, &ic, weight, INSERT_VALUES);
1205: }
1206: /* Elements */
1207: if (dof[3] == 1 && exf < startexf + nexf && eyf < starteyf + neyf && ezf < startezf + nezf) {
1208: DMStagStencil rowc, colf;
1209: PetscInt ir, ic;
1210: const PetscScalar weight = 0.125;
1212: colf.i = exf;
1213: colf.j = eyf;
1214: colf.k = ezf;
1215: colf.c = 0;
1216: colf.loc = DMSTAG_ELEMENT;
1217: rowc.i = exc;
1218: rowc.j = eyc;
1219: rowc.k = ezc;
1220: rowc.c = 0;
1221: rowc.loc = DMSTAG_ELEMENT;
1222: DMStagStencilToIndexLocal(dmc, dim, 1, &rowc, &ir);
1223: DMStagStencilToIndexLocal(dmf, dim, 1, &colf, &ic);
1224: MatSetValuesLocal(A, 1, &ir, 1, &ic, &weight, INSERT_VALUES);
1225: }
1226: }
1227: }
1228: }
1229: return 0;
1230: }