Actual source code: dmakkt.c
petsc-3.3-p0 2012-06-05
1: #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/
2: #include <../src/dm/impls/akkt/dmakkt.h> /*I "petscdmakkt.h" I*/
3: #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/
5: #undef __FUNCT__
7: PetscErrorCode DMAKKTSetDM(DM dm, DM ddm) {
8: PetscBool iskkt;
9: DM_AKKT *kkt = (DM_AKKT*)(dm->data);
14: PetscObjectTypeCompare((PetscObject)dm, DMAKKT, &iskkt);
15: if(!iskkt) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM not of type DMAKKT");
16: PetscObjectReference((PetscObject)ddm);
17: if(kkt->dm) {
18: DMDestroy(&(kkt->dm));
19: }
20: kkt->dm = ddm;
21: DMDestroy(&(kkt->cdm));
22: MatDestroy(&(kkt->Pfc));
23: dm->setupcalled = PETSC_FALSE;
24:
25: return(0);
26: }
28: #undef __FUNCT__
30: PetscErrorCode DMAKKTGetDM(DM dm, DM *ddm) {
31: PetscBool iskkt;
32: DM_AKKT *kkt = (DM_AKKT*)(dm->data);
36: PetscObjectTypeCompare((PetscObject)dm, DMAKKT, &iskkt);
37: if(!iskkt) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM not of type DMAKKT");
38: if(ddm) {
39: PetscObjectReference((PetscObject)(kkt->dm));
40: *ddm = kkt->dm;
41: }
42: return(0);
43: }
45: #undef __FUNCT__
47: PetscErrorCode DMAKKTSetMatrix(DM dm, Mat Aff) {
48: PetscBool iskkt;
50: DM_AKKT *kkt = (DM_AKKT*)(dm->data);
54: PetscObjectTypeCompare((PetscObject)dm, DMAKKT, &iskkt);
55: if(!iskkt) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM not of type DMAKKT");
56: PetscObjectReference((PetscObject)Aff);
57: if(kkt->Aff) {
58: MatDestroy(&(kkt->Aff));
59: }
60: kkt->Aff = Aff;
62: DMDestroy(&(kkt->cdm));
63: MatDestroy(&(kkt->Pfc));
64: dm->setupcalled = PETSC_FALSE;
65: return(0);
66: }
68: #undef __FUNCT__
70: PetscErrorCode DMAKKTGetMatrix(DM dm, Mat *Aff) {
71: PetscBool iskkt;
73: DM_AKKT *kkt = (DM_AKKT*)(dm->data);
76: PetscObjectTypeCompare((PetscObject)dm, DMAKKT, &iskkt);
77: if(!iskkt) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM not of type DMAKKT");
78: if(Aff) {
79: PetscObjectReference((PetscObject)(kkt->Aff));
80: *Aff = kkt->Aff;
81: }
82: return(0);
83: }
85: #undef __FUNCT__
87: PetscErrorCode DMAKKTSetFieldDecompositionName(DM dm, const char* dname) {
88: PetscBool iskkt;
89: DM_AKKT *kkt = (DM_AKKT*)(dm->data);
94: PetscObjectTypeCompare((PetscObject)dm, DMAKKT, &iskkt);
95: if(!iskkt) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM not of type DMAKKT");
96: if(kkt->dname) {
97: PetscStrncpy(kkt->dname, dname, DMAKKT_DECOMPOSITION_NAME_LEN);
98: }
99: DMDestroy(&(kkt->cdm));
100: MatDestroy(&(kkt->Pfc));
101: dm->setupcalled = PETSC_FALSE;
102:
103: return(0);
104: }
106: #undef __FUNCT__
108: PetscErrorCode DMAKKTGetFieldDecompositionName(DM dm, char** dname) {
109: PetscBool iskkt;
110: DM_AKKT *kkt = (DM_AKKT*)(dm->data);
115: PetscObjectTypeCompare((PetscObject)dm, DMAKKT, &iskkt);
116: if(!iskkt) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM not of type DMAKKT");
117: if(dname) {
118: *dname = PETSC_NULL;
119: if(kkt->dname) {
120: PetscStrallocpy(kkt->dname, dname);
121: }
122: }
123: return(0);
124: }
127: #undef __FUNCT__
129: PetscErrorCode DMAKKTSetFieldDecomposition(DM dm, PetscInt n, const char* const *names, IS *iss, DM *dms) {
130: PetscBool iskkt;
131: DM_AKKT *kkt = (DM_AKKT*)(dm->data);
132: PetscInt i;
139: PetscObjectTypeCompare((PetscObject)dm, DMAKKT, &iskkt);
140: if(!iskkt) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM not of type DMAKKT");
141: if(n < 1 || n > 2) SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Number of parts in decomposition must be between 1 and 2. Got %D instead",n);
142: for(i = 0; i < 2; ++i) {
143: if(kkt->names[i]) {
144: PetscFree(kkt->names[i]);
145: }
146: if(names[i]){
147: PetscStrallocpy(names[i], &(kkt->names[i]));
148: }
149: if(iss[i]) {
150: PetscObjectReference((PetscObject)iss[i]);
151: }
152: if(kkt->isf[i]) {
153: ISDestroy(&(kkt->isf[i]));
154: }
155: kkt->isf[i] = iss[i];
156: if(dms[i]) {
157: PetscObjectReference((PetscObject)dms[i]);
158: }
159: if(kkt->dmf[i]) {
160: DMDestroy(&(kkt->dmf[i]));
161: }
162: kkt->dmf[i] = dms[i];
163: }
164: DMDestroy(&(kkt->cdm));
165: MatDestroy(&(kkt->Pfc));
166: dm->setupcalled = PETSC_FALSE;
167:
168: return(0);
169: }
171: #undef __FUNCT__
173: PetscErrorCode DMAKKTGetFieldDecomposition(DM dm, PetscInt *n, char*** names, IS **iss, DM **dms) {
174: PetscBool iskkt;
175: DM_AKKT *kkt = (DM_AKKT*)(dm->data);
176: PetscInt i;
183: PetscObjectTypeCompare((PetscObject)dm, DMAKKT, &iskkt);
184: if(!iskkt) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM not of type DMAKKT");
185: if(n) *n = 2;
186: if(names) {
187: if(kkt->names) {
188: PetscMalloc(sizeof(char*)*2, names);
189: }
190: else {
191: *names = PETSC_NULL;
192: }
193: }
194: if(iss) {
195: if(kkt->isf) {
196: PetscMalloc(sizeof(IS)*2, iss);
197: }
198: else {
199: *iss = PETSC_NULL;
200: }
201: }
202: if(dms) {
203: if(kkt->dmf) {
204: PetscMalloc(sizeof(DM)*2, dms);
205: }
206: else {
207: *dms = PETSC_NULL;
208: }
209: }
210: for(i = 0; i < 2; ++i) {
211: if(names && kkt->names){
212: PetscStrallocpy(kkt->names[i],(*names)+i);
213: }
214: if(iss && kkt->isf) {
215: PetscObjectReference((PetscObject)kkt->isf[i]);
216: (*iss)[i] = kkt->isf[i];
217: }
218: if(dms && kkt->dmf) {
219: PetscObjectReference((PetscObject)kkt->dmf[i]);
220: (*dms)[i] = kkt->dmf[i];
221: }
222: }
223: return(0);
224: }
227: #undef __FUNCT__
229: PetscErrorCode DMSetFromOptions_AKKT(DM dm) {
231: DM_AKKT* kkt = (DM_AKKT*)(dm->data);
233: PetscOptionsBool("-dm_akkt_duplicate_mat",
234: "Duplicate underlying Mat in DMCreateMatrix",
235: "DMAKKTSetDupulicateMat",
236: kkt->duplicate_mat,
237: &kkt->duplicate_mat,
238: PETSC_NULL);
239:
240: PetscOptionsBool("-dm_akkt_detect_saddle_point",
241: "Identify dual variables by zero diagonal entries",
242: "DMAKKTSetDetectSaddlePoint",
243: kkt->detect_saddle_point,
244: &kkt->detect_saddle_point,
245: PETSC_NULL);
246:
247: PetscOptionsString("-dm_akkt_decomposition_name",
248: "Name of primal-dual decomposition to request from DM",
249: "DMAKKTSetFieldDecompositionName",
250: kkt->dname,
251: kkt->dname,
252: DMAKKT_DECOMPOSITION_NAME_LEN,
253: PETSC_NULL);
254:
255: dm->setupcalled = PETSC_FALSE;
256: return(0);
257: }
259: #undef __FUNCT__
261: PetscErrorCode DMSetUp_AKKT(DM dm) {
262: DM_AKKT *kkt = (DM_AKKT*)(dm->data);
265: if(dm->setupcalled) return(0);
266: if(!kkt->Aff){
267: if(kkt->dm) {
268: DMCreateMatrix(kkt->dm, MATAIJ, &kkt->Aff);
269: }
270: else SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Neither matrix nor DM set");
271: }
272: if(!kkt->isf[0] && !kkt->isf[0]) {
273: if(kkt->detect_saddle_point) {
274: MatFindZeroDiagonals(kkt->Aff,&kkt->isf[1]);
275: }
276: else if(kkt->dm && kkt->dname) {
277: DM ddm;
278: PetscInt n;
279: char **names;
280: IS *iss;
281: DM *dms;
282: PetscInt i;
283: DMCreateFieldDecompositionDM(kkt->dm, kkt->dname, &ddm);
284: DMCreateFieldDecomposition(ddm, &n, &names, &iss, &dms);
285: if(n < 1 || n > 2)
286: SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Number of parts in decomposition %s must be between 1 and 2. Got %D instead",kkt->dname, n);
287: for(i = 0; i < n; ++i) {
288: if(!iss[i] && dms[i]) {
289: const char* label;
290: if(i == 0)
291: label = "primal";
292: else
293: label = "dual";
294: SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Decomposition defines %s subDM, but no embedding IS is given", label);
295: }
296: }
297: DMAKKTSetFieldDecomposition(dm, n, (const char**)names, iss, dms);
298: for(i = 0; i < n; ++i) {
299: PetscFree(names[i]);
300: ISDestroy(&(iss[i]));
301: DMDestroy(&(dms[i]));
302: }
303: }
304: }
305: if(!kkt->isf[0] && !kkt->isf[1]) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Decomposition ISs not set and could not be derived. ");
306: if(!kkt->isf[0] || !kkt->isf[1]) {
307: PetscInt lstart, lend;
308: MatGetOwnershipRange(kkt->Aff, &lstart, &lend);
309: if(!kkt->isf[0]) {
310: ISComplement(kkt->isf[0], lstart, lend, kkt->isf+1);
311: }
312: else {
313: ISComplement(kkt->isf[1], lstart, lend, kkt->isf+0);
314: }
315: }
316: /* FIX: Should we allow a combination of empty kkt->dmf[0] and non-empty kkt->dmf[1]? */
317: if(!kkt->dmf[0]) {
318: /* Construct a GAMG proxy to coarsen the primal block. */
319: Mat A0f0f;
320: IS is00;
321: PetscInt lstart, lend;
322: const char* primal = {"all"};
323: DMCreate(((PetscObject)dm)->comm, kkt->dmf+0);
324: DMSetType(kkt->dmf[0],DMAKKT);
325: MatGetSubMatrix(kkt->Aff, kkt->isf[0], kkt->isf[0], MAT_INITIAL_MATRIX, &A0f0f);
326: DMAKKTSetMatrix(kkt->dmf[0], A0f0f);
327: MatGetOwnershipRange(A0f0f, &lstart, &lend);
328: ISCreateStride(((PetscObject)A0f0f)->comm, lend-lstart, lstart, 1, &is00);
329: DMAKKTSetFieldDecomposition(kkt->dmf[0], 1, &primal, &is00, PETSC_NULL);
330: }
331: dm->setupcalled = PETSC_TRUE;
332: return(0);
333: }
335: /*
336: This routine will coarsen the 11-block only using the 00-block prolongation (P0f0c), the 10 block and GAMG.
337: The result is the 11-block prolongator (P1f1c).
338: */
339: #undef __FUNCT__
341: PetscErrorCode DMCoarsen_AKKT_GAMG11(DM dm, Mat P0f0c, Mat *P1f1c_out) {
343: DM_AKKT* kkt = (DM_AKKT*)(dm->data);
344: Mat Aff = kkt->Aff; /* fine-level KKT matrix */
345: Mat A1f0f; /* fine-level dual (constraint) Jacobian */
346: Mat A1f0c; /* = A1f0f*P0f0c coarsen only primal indices */
347: Mat B1f1f; /* = A1f0c'*A1f0c */
348: PC gamg11;/* Use PCGAMG internally to get access to some of its methods to operate on B1f1f = A1f0c*A1f0c', where A1f0c = A1f0f*P0f0c. */
349: PC_GAMG* pc_gamg11;
350: Mat G1f1f; /* = Graph(B1f1f) */
351: Mat P1f1c; /* = Prolongator(G1f1f); */
352: PetscCoarsenData *coarsening;
355: /*
356: What is faster:
357: - A0c1f = P0f0c'*A0f1f followed by B1f1f = A0c1f'*A0c1f, or
358: - A1f0c = A1f0f*P0f0c followed by B1f1f = A1f0c*A1f0c'?
359: My bet is on the latter:
360: - fewer transpositions inside MatMatMult and row indices are always local.
361: */
363: MatGetSubMatrix(Aff, kkt->isf[1], kkt->isf[0], MAT_INITIAL_MATRIX, &A1f0f);
364: if(kkt->transposeP) {
365: MatMatTransposeMult(A1f0f,P0f0c,MAT_INITIAL_MATRIX, PETSC_DEFAULT, &A1f0c);
366: }
367: MatMatMult(A1f0f,P0f0c,MAT_INITIAL_MATRIX, PETSC_DEFAULT, &A1f0c);
368: MatMatTransposeMult(A1f0c, A1f0c, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &B1f1f);
369:
370: /* We create PCGAMG here since it is only needed for coarsening and we don't want to have to carry the attendant data structures, if we don't need them. */
371: PCCreate(((PetscObject)dm)->comm, &gamg11);
372: /* This must be an aggregating GAMG. */
373: PCSetType(gamg11, PCGAMG);
374: PCGAMGSetSquareGraph(gamg11, PETSC_FALSE);
375: /*
376: Observe that we want to "square" A1f0c before passing it (B1f1f) to GAMG.
377: This is not because we are not sure how GAMG will deal with a (potentially) non-square matrix,
378: but rather because if we asked GAMG to square it, it would also smooth the resulting prolongator.
379: At least PC_GAMG_AGG would, and we need an unsmoothed prolongator.
380: */
381: PCSetOperators(gamg11, B1f1f, B1f1f, DIFFERENT_NONZERO_PATTERN);
382: /* FIX: Currently there is no way to tell GAMG to coarsen onto a give comm, but it shouldn't be hard to hack that stuff in. */
383: pc_gamg11 = (PC_GAMG*)(gamg11->data);
384: pc_gamg11->graph(gamg11, B1f1f, &G1f1f);
385: pc_gamg11->coarsen(gamg11, &G1f1f, &coarsening);
386: pc_gamg11->prolongator(gamg11, B1f1f, G1f1f, coarsening, &P1f1c);
388: MatDestroy(&A1f0f);
389: MatDestroy(&A1f0c);
390: MatDestroy(&B1f1f);
391: MatDestroy(&G1f1f);
392: PCDestroy(&gamg11);
393:
394: *P1f1c_out = P1f1c;
396: return(0);
397: }
399: #undef __FUNCT__
401: PetscErrorCode DMCoarsen_AKKT(DM dm, MPI_Comm comm, DM *cdm) {
403: DM_AKKT* kkt = (DM_AKKT*)(dm->data);
404: Mat Acc; /* coarse-level KKT matrix */
405: Mat P0f0c, P1f1c; /* Primal and dual block prolongator */
406: DM dmc[2] = {PETSC_NULL, PETSC_NULL}; /* Coarse subDMs defining the block prolongators and the coarsened decomposition. */
407: PetscInt M0,N0,M1,N1; /* Sizes of P0f0c and P1f1c. */
408: PetscInt start0,end0,start1,end1; /* Ownership ranges for P0f0c and P1f1c. */
409: static Mat mats[4] = {PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL}; /* Used to construct MatNest out of pieces. */
410: IS isc[2]; /* Used to construct MatNest out of pieces and to define the coarsened decomposition. */
412: if(!cdm) return(0);
413: if(kkt->cdm) {
414: PetscObjectReference((PetscObject)(kkt->cdm));
415: *cdm = kkt->cdm;
416: return(0);
417: }
418: /* Coarsen the 00 block with the attached DM and obtain the primal prolongator. */
419: if(kkt->dmf[0]) {
420: DMCoarsen(kkt->dmf[0], comm, dmc+0);
421: DMCreateInterpolation(dmc[0], kkt->dmf[0], &P0f0c, PETSC_NULL);
422: }
423: else SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Could not coarsen the primal block: primal subDM not set.");
424:
425: /* Should P0f0c be transposed to act as a prolongator (i.e., to map from coarse to fine). */
426: MatGetSize(P0f0c, &M0, &N0);
427: if(M0 == N0) SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE,"Primal prolongator is square with size %D: cannot distinguish coarse from fine",M0);
428: if(M0 < N0) kkt->transposeP = PETSC_TRUE;
429: else kkt->transposeP = PETSC_FALSE;
430: /* See if the 11 block can be coarsened with an attached DM. If so, we are done. Otherwise, use GAMG to coarsen 11. */
431: if(kkt->dmf[1]) {
432: DMCoarsen(kkt->dmf[1], comm, &dmc[1]);
433: DMCreateInterpolation(dmc[1], kkt->dmf[1], &P1f1c, PETSC_NULL);
434: }
435: else {
436: DMCoarsen_AKKT_GAMG11(dm, P0f0c, &P1f1c);
437: }
438: /* Determine whether P1f1c should be transposed in order to act as a prolongator (i.e., to map from coarse to fine). */
439: MatGetSize(P1f1c, &M1, &N1);
440: if(M1 == N1) SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Dual prlongator is square with size %D: cannot distinguish coarse from fine", M1);
441: if((M1 < N1 && !kkt->transposeP) || (M1 >= N1 && kkt->transposeP)) {
442: Mat P1f1ct;
443: MatTranspose(P1f1c, MAT_INITIAL_MATRIX, &P1f1ct);
444: MatDestroy(&P1f1c);
445: P1f1c = P1f1ct;
446: }
447: /* MatNest P0f0c, P1f1c together into Pfc. */
448: mats[0] = P0f0c; mats[3] = P1f1c;
449: MatGetOwnershipRange(P0f0c, &start0, &end0);
450: MatGetOwnershipRange(P1f1c, &start1, &end1);
451: ISCreateStride(((PetscObject)dm)->comm, end0-start0,start0,1,isc+0);
452: ISCreateStride(((PetscObject)dm)->comm, end1-start1,start1,1,isc+1);
453: if(kkt->transposeP) {
454: MatCreateNest(((PetscObject)dm)->comm,2,isc,2,kkt->isf,mats,&(kkt->Pfc));
455: }
456: else {
457: MatCreateNest(((PetscObject)dm)->comm,2,kkt->isf,2,isc,mats,&(kkt->Pfc));
458: }
459: MatDestroy(&P0f0c);
460: MatDestroy(&P1f1c);
461: /* Coarsening the underlying matrix and primal-dual decomposition. */
462: /*
463: We do not coarsen the underlying DM because
464: (a) Its coarsening may be incompatible with the specialized KKT-aware coarsening of the blocks defined here.
465: (b) Even if the coarsening the decomposition is compatible with the decomposition of the coarsening, we can
466: pick the former without loss of generality.
467: (c) Even if (b) is true, the embeddings (IS) of the coarsened subDMs are potentially different now from what
468: they would be in the coarsened DM; thus, embeddings would have to be supplied manually anyhow.
469: (d) In the typical situation we should only use the primal subDM for coarsening -- the whole point of
470: DMAKKT is that the dual block coarsening should be derived from the primal block coarsening for compatibility.
471: If we are given both subDMs, DMAKKT essentially becomes a version of DMComposite, in which case the composition
472: of the coarsened decomposition is by definition the coarsening of the whole system DM.
473: */
474: /* Create the coarser DM. */
475: DMCreate(((PetscObject)dm)->comm, &(kkt->cdm));
476: DMSetType(kkt->cdm, DMAKKT);
477: /* Coarsen the underlying matrix. */
478: MatPtAP(kkt->Aff, kkt->Pfc, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Acc);
479: DMAKKTSetMatrix(dm, Acc);
480: /* Set the coarsened decomposition. */
481: DMAKKTSetFieldDecomposition(kkt->cdm, 2, (const char**)kkt->names, isc, dmc);
482: return(0);
483: }
485: #undef __FUNCT__
487: PetscErrorCode DMCreateInterpolation_AKKT(DM cdm, DM fdm, Mat *interp, Vec* rscale) {
488: PetscBool iskkt;
490: DM_AKKT* fkkt = (DM_AKKT*)(fdm->data);
494: PetscObjectTypeCompare((PetscObject)cdm, DMAKKT, &iskkt);
495: if(!iskkt) SETERRQ(((PetscObject)cdm)->comm, PETSC_ERR_ARG_WRONG, "Coarse DM not of type DMAKKT");
496: PetscObjectTypeCompare((PetscObject)fdm, DMAKKT, &iskkt);
497: if(!iskkt) SETERRQ(((PetscObject)fdm)->comm, PETSC_ERR_ARG_WRONG, "Fine DM not of type DMAKKT");
498: if(fkkt->cdm != cdm) SETERRQ(((PetscObject)cdm)->comm, PETSC_ERR_ARG_WRONG, "Coarse DM must be obtained from fine via DMCoarsen");
499: if(interp) {
500: PetscObjectReference((PetscObject)(fkkt->Pfc));
501: *interp = fkkt->Pfc;
502: }
503: if(rscale) *rscale = PETSC_NULL;
505: return(0);
506: }
508: #undef __FUNCT__
510: PetscErrorCode DMCreateMatrix_AKKT(DM dm, const MatType type, Mat *A) {
511: PetscBool iskkt;
513: DM_AKKT* kkt = (DM_AKKT*)(dm->data);
516: PetscObjectTypeCompare((PetscObject)dm, DMAKKT, &iskkt);
517: if(!iskkt) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM not of type DMAKKT");
518: if(A) {
519: if(kkt->duplicate_mat) {
520: MatDuplicate(kkt->Aff, MAT_SHARE_NONZERO_PATTERN, A);
521: }
522: else {
523: PetscObjectReference((PetscObject)(kkt->Aff));
524: *A = kkt->Aff;
525: }
526: }
527: return(0);
528: }
530: #undef __FUNCT__
532: PetscErrorCode DMCreateGlobalVector_AKKT(DM dm, Vec *v) {
533: PetscBool iskkt;
535: DM_AKKT* kkt = (DM_AKKT*)(dm->data);
538: PetscObjectTypeCompare((PetscObject)dm, DMAKKT, &iskkt);
539: if(!iskkt) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM not of type DMAKKT");
540: if(v) {
541: MatGetVecs(kkt->Aff, v, PETSC_NULL);
542: }
543: return(0);
544: }
546: #undef __FUNCT__
548: PetscErrorCode DMDestroy_AKKT(DM dm) {
549: DM_AKKT *kkt = (DM_AKKT*)(dm->data);
550: PetscInt i;
553: MatDestroy(&(kkt->Aff));
554: DMDestroy(&(kkt->dm));
555: for(i = 0; i < 2; ++i) {
556: DMDestroy(&(kkt->dmf[i]));
557: ISDestroy(&(kkt->isf[i]));
558: PetscFree(kkt->names[i]);
559: }
560: DMDestroy(&(kkt->cdm));
561: MatDestroy(&(kkt->Pfc));
562: PetscFree(kkt);
563: return(0);
564: }
566: #undef __FUNCT__
568: PetscErrorCode DMView_AKKT(DM dm, PetscViewer v) {
569: DM_AKKT* kkt = (DM_AKKT*)(dm->data);
571: PetscBool isascii;
572: PetscInt i, tab, vtab;
573: const char* name, *prefix;
575: PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii);
576: if(!isascii) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_SUP, "No support for non-ASCII viewers");
577: PetscObjectGetTabLevel((PetscObject)dm, &tab);
578: PetscObjectGetName((PetscObject)dm, &name);
579: PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix);
580: PetscViewerASCIIUseTabs(v,PETSC_TRUE);
581: PetscViewerASCIIGetTab(v,&vtab);
582: PetscViewerASCIISetTab(v,tab);
583: PetscViewerASCIIPrintf(v, "DM Algebraic KKT, name: %s, prefix: %s\n", ((PetscObject)dm)->name, ((PetscObject)dm)->prefix);
584: if(kkt->dm) {
585: PetscViewerASCIIPrintf(v, "DM:\n");
586: PetscViewerASCIIPushTab(v);
587: DMView(kkt->dm,v);
588: PetscViewerASCIIPopTab(v);
589: }
590: if(kkt->Aff) {
591: PetscViewerASCIIPrintf(v, "Aff:\n");
592: PetscViewerASCIIPushTab(v);
593: MatView(kkt->Aff,v);
594: PetscViewerASCIIPopTab(v);
595: }
596: if(kkt->dname) {
597: PetscViewerASCIIPrintf(v, "Decomposition, name %s:\n");
598: }
599: for(i = 0; i < 2; ++i) {
600: const char* label;
601: if(i == 0) {
602: label = "Primal";
603: }
604: else {
605: label = "Dual";
606: }
607: if(kkt->names[i]) {
608: PetscViewerASCIIPrintf(v, "%s, name %s:\n", label, kkt->names[i]);
609: }
610: if(kkt->isf[i]){
611: PetscViewerASCIIPrintf(v, "%s, IS:\n",label);
612: PetscViewerASCIIPushTab(v);
613: ISView(kkt->isf[i],v);
614: PetscViewerASCIIPopTab(v);
615: }
616: if(kkt->dmf[i]){
617: PetscViewerASCIIPrintf(v, "%s, DM:\n", label);
618: PetscViewerASCIIPushTab(v);
619: DMView(kkt->dmf[i],v);
620: PetscViewerASCIIPopTab(v);
621: }
622: }
623: if(kkt->Pfc) {
624: PetscViewerASCIIPrintf(v, "Prolongation:\n");
625: PetscViewerASCIIPushTab(v);
626: MatView(kkt->Pfc,v);
627: PetscViewerASCIIPopTab(v);
628: }
629:
630: PetscViewerASCIISetTab(v,vtab);
631:
632: return(0);
633: }
635: #undef __FUNCT__
637: PetscErrorCode DMCreate_AKKT(DM dm) {
638: DM_AKKT *kkt;
641: PetscNewLog(dm, DM_AKKT, &kkt);
642: dm->data = kkt;
643: PetscObjectChangeTypeName((PetscObject)kkt,DMAKKT);
644: kkt->dmf[0] = kkt->dmf[1] = PETSC_NULL;
645: kkt->isf[0] = kkt->isf[1] = PETSC_NULL;
646: kkt->names[0] = kkt->names[1] = PETSC_NULL;
648: kkt->duplicate_mat = PETSC_FALSE;
649: kkt->detect_saddle_point = PETSC_FALSE;
651: dm->ops->createglobalvector = DMCreateGlobalVector_AKKT;
652: dm->ops->creatematrix = DMCreateMatrix_AKKT;
653: dm->ops->createinterpolation = DMCreateInterpolation_AKKT;
654: dm->ops->coarsen = DMCoarsen_AKKT;
655: dm->ops->destroy = DMDestroy_AKKT;
656: dm->ops->view = DMView_AKKT;
657: dm->ops->setfromoptions = DMSetFromOptions_AKKT;
658: dm->ops->setup = DMSetUp_AKKT;
660: return(0);
661: }