Actual source code: forest.c
petsc-3.7.3 2016-08-01
1: #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/
2: #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/
3: #include <petscsf.h>
5: PetscBool DMForestPackageInitialized = PETSC_FALSE;
7: typedef struct _DMForestTypeLink *DMForestTypeLink;
9: struct _DMForestTypeLink
10: {
11: char *name;
12: DMForestTypeLink next;
13: };
15: DMForestTypeLink DMForestTypeList;
19: static PetscErrorCode DMForestPackageFinalize(void)
20: {
21: DMForestTypeLink oldLink, link = DMForestTypeList;
25: while (link) {
26: oldLink = link;
27: PetscFree(oldLink->name);
28: link = oldLink->next;
29: PetscFree(oldLink);
30: }
31: return(0);
32: }
36: static PetscErrorCode DMForestPackageInitialize(void)
37: {
41: if (DMForestPackageInitialized) return(0);
42: DMForestPackageInitialized = PETSC_TRUE;
43: DMForestRegisterType(DMFOREST);
44: PetscRegisterFinalize(DMForestPackageFinalize);
45: return(0);
46: }
50: /*@C
51: DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct)
53: Not Collective
55: Input parameter:
56: . name - the name of the type
58: Level: advanced
60: .seealso: DMFOREST, DMIsForest()
61: @*/
62: PetscErrorCode DMForestRegisterType(DMType name)
63: {
64: DMForestTypeLink link;
68: DMForestPackageInitialize();
69: PetscNew(&link);
70: PetscStrallocpy(name,&link->name);
71: link->next = DMForestTypeList;
72: DMForestTypeList = link;
73: return(0);
74: }
78: /*@
79: DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes
81: Not Collective
83: Input parameter:
84: . dm - the DM object
86: Output parameter:
87: . isForest - whether dm is a subtype of DMFOREST
89: Level: intermediate
91: .seealso: DMFOREST, DMForestRegisterType()
92: @*/
93: PetscErrorCode DMIsForest(DM dm, PetscBool *isForest)
94: {
95: DMForestTypeLink link = DMForestTypeList;
99: while (link) {
100: PetscBool sameType;
101: PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);
102: if (sameType) {
103: *isForest = PETSC_TRUE;
104: return(0);
105: }
106: link = link->next;
107: }
108: *isForest = PETSC_FALSE;
109: return(0);
110: }
114: /*@
115: DMForestTemplate - Create a new DM that will be adapted from a source DM. The new DM reproduces the configuration
116: of the source, but is not yet setup, so that the user can then define only the ways that the new DM should differ
117: (by, e.g., refinement or repartitioning). The source DM is also set as the adaptivity source DM of the new DM (see
118: DMForestSetAdaptivityForest()).
120: Collective on dm
122: Input Parameters:
123: + dm - the source DM object
124: - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen())
126: Output Parameter:
127: . tdm - the new DM object
129: Level: intermediate
131: .seealso: DMForestSetAdaptivityForest()
132: @*/
133: PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm)
134: {
135: DM_Forest *forest = (DM_Forest *) dm->data;
136: DMType type;
137: DM base;
138: DMForestTopology topology;
139: PetscInt dim, overlap, ref, factor;
140: DMForestAdaptivityStrategy strat;
141: PetscDS ds;
142: void *ctx;
143: PetscErrorCode ierr;
147: DMCreate(PetscObjectComm((PetscObject)dm),tdm);
148: DMGetType(dm,&type);
149: DMSetType(*tdm,type);
150: DMForestGetBaseDM(dm,&base);
151: DMForestSetBaseDM(*tdm,base);
152: DMForestGetTopology(dm,&topology);
153: DMForestSetTopology(*tdm,topology);
154: DMForestGetAdjacencyDimension(dm,&dim);
155: DMForestSetAdjacencyDimension(*tdm,dim);
156: DMForestGetPartitionOverlap(dm,&overlap);
157: DMForestSetPartitionOverlap(*tdm,overlap);
158: DMForestGetMinimumRefinement(dm,&ref);
159: DMForestSetMinimumRefinement(*tdm,ref);
160: DMForestGetMaximumRefinement(dm,&ref);
161: DMForestSetMaximumRefinement(*tdm,ref);
162: DMForestGetAdaptivityStrategy(dm,&strat);
163: DMForestSetAdaptivityStrategy(*tdm,strat);
164: DMForestGetGradeFactor(dm,&factor);
165: DMForestSetGradeFactor(*tdm,factor);
166: if (forest->ftemplate) {
167: (forest->ftemplate) (dm, *tdm);
168: }
169: DMForestSetAdaptivityForest(*tdm,dm);
170: DMGetDS(dm,&ds);
171: DMSetDS(*tdm,ds);
172: DMGetApplicationContext(dm,&ctx);
173: DMSetApplicationContext(*tdm,&ctx);
174: if (dm->maxCell) {
175: const PetscReal *maxCell, *L;
176: const DMBoundaryType *bd;
178: DMGetPeriodicity(dm,&maxCell,&L,&bd);
179: DMSetPeriodicity(*tdm,maxCell,L,bd);
180: }
181: return(0);
182: }
184: static PetscErrorCode DMInitialize_Forest(DM dm);
188: PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm)
189: {
190: DM_Forest *forest = (DM_Forest *) dm->data;
191: const char *type;
195: forest->refct++;
196: (*newdm)->data = forest;
197: PetscObjectGetType((PetscObject) dm, &type);
198: PetscObjectChangeTypeName((PetscObject) *newdm, type);
199: DMInitialize_Forest(*newdm);
200: return(0);
201: }
205: static PetscErrorCode DMDestroy_Forest(DM dm)
206: {
207: DM_Forest *forest = (DM_Forest*) dm->data;
211: if (--forest->refct > 0) return(0);
212: if (forest->destroy) {forest->destroy(dm);}
213: PetscSFDestroy(&forest->cellSF);
214: PetscSFDestroy(&forest->preCoarseToFine);
215: PetscSFDestroy(&forest->coarseToPreFine);
216: PetscFree(forest->adaptLabel);
217: PetscFree(forest->adaptStrategy);
218: DMDestroy(&forest->base);
219: DMDestroy(&forest->adapt);
220: PetscFree(forest->topology);
221: PetscFree(forest);
222: return(0);
223: }
227: /*@C
228: DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase. The topology is a string (e.g.
229: "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest durint
230: DMSetUp().
232: Logically collective on dm
234: Input parameters:
235: + dm - the forest
236: - topology - the topology of the forest
238: Level: intermediate
240: .seealso(): DMForestGetTopology(), DMForestSetBaseDM()
241: @*/
242: PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
243: {
244: DM_Forest *forest = (DM_Forest *) dm->data;
249: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup");
250: PetscFree(forest->topology);
251: PetscStrallocpy((const char *)topology,(char **) &forest->topology);
252: return(0);
253: }
257: /*@C
258: DMForestGetTopology - Get a string describing the topology of a DMForest.
260: Not collective
262: Input parameter:
263: . dm - the forest
265: Output parameter:
266: . topology - the topology of the forest (e.g., 'cube', 'shell')
268: Level: intermediate
270: .seealso: DMForestSetTopology()
271: @*/
272: PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
273: {
274: DM_Forest *forest = (DM_Forest *) dm->data;
279: *topology = forest->topology;
280: return(0);
281: }
285: /*@
286: DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest. The
287: forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its
288: base. In general, two forest must share a bse to be comparable, to do things like construct interpolators.
290: Logically collective on dm
292: Input Parameters:
293: + dm - the forest
294: - base - the base DM of the forest
296: Level: intermediate
298: .seealso(): DMForestGetBaseDM()
299: @*/
300: PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
301: {
302: DM_Forest *forest = (DM_Forest *) dm->data;
303: PetscInt dim, dimEmbed;
308: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup");
309: PetscObjectReference((PetscObject)base);
310: DMDestroy(&forest->base);
311: forest->base = base;
312: if (base) {
314: DMGetDimension(base,&dim);
315: DMSetDimension(dm,dim);
316: DMGetCoordinateDim(base,&dimEmbed);
317: DMSetCoordinateDim(dm,dimEmbed);
318: }
319: return(0);
320: }
324: /*@
325: DMForestGetBaseDM - Get the base DM of a DMForest forest. The forest will be hierarchically refined from the base,
326: and all refinements/coarsenings of the forest will share its base. In general, two forest must share a bse to be
327: comparable, to do things like construct interpolators.
329: Not collective
331: Input Parameter:
332: . dm - the forest
334: Output Parameter:
335: . base - the base DM of the forest
337: Level: intermediate
339: .seealso(); DMForestSetBaseDM()
340: @*/
341: PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
342: {
343: DM_Forest *forest = (DM_Forest *) dm->data;
348: *base = forest->base;
349: return(0);
350: }
354: /*@
355: DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
356: adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp(). Usually not needed
357: by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this
358: routine.
360: Logically collective on dm
362: Input Parameter:
363: + dm - the new forest, which will be constructed from adapt
364: - adapt - the old forest
366: Level: intermediate
368: .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose()
369: @*/
370: PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt)
371: {
372: DM_Forest *forest;
373: PetscErrorCode ierr;
378: forest = (DM_Forest *) dm->data;
379: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
380: switch (forest->adaptPurpose) {
381: case DM_FOREST_KEEP:
382: PetscObjectReference((PetscObject)adapt);
383: DMDestroy(&(forest->adapt));
384: forest->adapt = adapt;
385: break;
386: case DM_FOREST_REFINE:
387: DMSetCoarseDM(dm,adapt);
388: break;
389: case DM_FOREST_COARSEN:
390: DMSetFineDM(dm,adapt);
391: break;
392: default:
393: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
394: }
395: return(0);
396: }
400: /*@
401: DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
403: Not collective
405: Input Parameter:
406: . dm - the forest
408: Output Parameter:
409: . adapt - the forest from which dm is/was adapted
411: Level: intermediate
413: .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
414: @*/
415: PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
416: {
417: DM_Forest *forest;
418: PetscErrorCode ierr;
422: forest = (DM_Forest *) dm->data;
423: switch (forest->adaptPurpose) {
424: case DM_FOREST_KEEP:
425: *adapt = forest->adapt;
426: break;
427: case DM_FOREST_REFINE:
428: DMGetCoarseDM(dm,adapt);
429: break;
430: case DM_FOREST_COARSEN:
431: DMGetFineDM(dm,adapt);
432: break;
433: default:
434: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
435: }
436: return(0);
437: }
441: /*@
442: DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
443: source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening
444: (DM_FOREST_COARSEN), or undefined (DM_FOREST_NONE). This only matters for the purposes of reference counting:
445: during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
446: relationship (see DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does
447: not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
448: cause a memory leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
450: Logically collective on dm
452: Input Parameters:
453: + dm - the forest
454: - purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN)
456: Level: advanced
458: .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
459: @*/
460: PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose purpose)
461: {
462: DM_Forest *forest;
466: forest = (DM_Forest *) dm->data;
467: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
468: if (purpose != forest->adaptPurpose) {
469: DM adapt;
471: DMForestGetAdaptivityForest(dm,&adapt);
472: PetscObjectReference((PetscObject)adapt);
473: DMForestSetAdaptivityForest(dm,NULL);
474: forest->adaptPurpose = purpose;
475: DMForestSetAdaptivityForest(dm,adapt);
476: DMDestroy(&adapt);
477: }
478: return(0);
479: }
483: /*@
484: DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
485: DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening (DM_FOREST_COARSEN), or
486: undefined (DM_FOREST_NONE). This only matters for the purposes of reference counting: during DMDestroy(), cyclic
487: references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
488: DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does not maintain a
489: reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
490: leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
492: Not collective
494: Input Parameter:
495: . dm - the forest
497: Output Parameter:
498: . purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN)
500: Level: advanced
502: .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
503: @*/
504: PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose *purpose)
505: {
506: DM_Forest *forest;
509: forest = (DM_Forest *) dm->data;
510: *purpose = forest->adaptPurpose;
511: return(0);
512: }
516: /*@
517: DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
518: cell adjacency (for the purposes of partitioning and overlap).
520: Logically collective on dm
522: Input Parameters:
523: + dm - the forest
524: - adjDim - default 0 (i.e., vertices determine adjacency)
526: Level: intermediate
528: .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap()
529: @*/
530: PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
531: {
532: PetscInt dim;
533: DM_Forest *forest = (DM_Forest *) dm->data;
534: PetscErrorCode ierr;
538: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
539: if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim);
540: DMGetDimension(dm,&dim);
541: if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim);
542: forest->adjDim = adjDim;
543: return(0);
544: }
548: /*@
549: DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
550: e.g., adjacency based on facets can be specified by codimension 1 in all cases)
552: Logically collective on dm
554: Input Parameters:
555: + dm - the forest
556: - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
558: Level: intermediate
560: .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
561: @*/
562: PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
563: {
564: PetscInt dim;
565: PetscErrorCode ierr;
569: DMGetDimension(dm,&dim);
570: DMForestSetAdjacencyDimension(dm,dim-adjCodim);
571: return(0);
572: }
576: /*@
577: DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
578: purposes of partitioning and overlap).
580: Not collective
582: Input Parameter:
583: . dm - the forest
585: Output Parameter:
586: . adjDim - default 0 (i.e., vertices determine adjacency)
588: Level: intermediate
590: .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap()
591: @*/
592: PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
593: {
594: DM_Forest *forest = (DM_Forest *) dm->data;
599: *adjDim = forest->adjDim;
600: return(0);
601: }
605: /*@
606: DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
607: e.g., adjacency based on facets can be specified by codimension 1 in all cases)
609: Not collective
611: Input Parameter:
612: . dm - the forest
614: Output Parameter:
615: . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
617: Level: intermediate
619: .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension()
620: @*/
621: PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
622: {
623: DM_Forest *forest = (DM_Forest *) dm->data;
624: PetscInt dim;
630: DMGetDimension(dm,&dim);
631: *adjCodim = dim - forest->adjDim;
632: return(0);
633: }
637: /*@
638: DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
639: partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
640: adjacent cells
642: Logically collective on dm
644: Input Parameters:
645: + dm - the forest
646: - overlap - default 0
648: Level: intermediate
650: .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
651: @*/
652: PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
653: {
654: DM_Forest *forest = (DM_Forest *) dm->data;
658: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
659: if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap);
660: forest->overlap = overlap;
661: return(0);
662: }
666: /*@
667: DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
668: > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
670: Not collective
672: Input Parameter:
673: . dm - the forest
675: Output Parameter:
676: . overlap - default 0
678: Level: intermediate
680: .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
681: @*/
682: PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
683: {
684: DM_Forest *forest = (DM_Forest *) dm->data;
689: *overlap = forest->overlap;
690: return(0);
691: }
695: /*@
696: DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
697: DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest
698: (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.
700: Logically collective on dm
702: Input Parameters:
703: + dm - the forest
704: - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
706: Level: intermediate
708: .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
709: @*/
710: PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
711: {
712: DM_Forest *forest = (DM_Forest *) dm->data;
716: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
717: forest->minRefinement = minRefinement;
718: return(0);
719: }
723: /*@
724: DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
725: DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest (see
726: DMForestGetAdaptivityForest()), this limits the amount of coarsening.
728: Not collective
730: Input Parameter:
731: . dm - the forest
733: Output Parameter:
734: . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
736: Level: intermediate
738: .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
739: @*/
740: PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
741: {
742: DM_Forest *forest = (DM_Forest *) dm->data;
747: *minRefinement = forest->minRefinement;
748: return(0);
749: }
753: /*@
754: DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
755: DM, see DMForestGetBaseDM()) allowed in the forest.
757: Logically collective on dm
759: Input Parameters:
760: + dm - the forest
761: - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
763: Level: intermediate
765: .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
766: @*/
767: PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
768: {
769: DM_Forest *forest = (DM_Forest *) dm->data;
773: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup");
774: forest->initRefinement = initRefinement;
775: return(0);
776: }
780: /*@
781: DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
782: DMForestGetBaseDM()) allowed in the forest.
784: Not collective
786: Input Parameter:
787: . dm - the forest
789: Output Paramater:
790: . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
792: Level: intermediate
794: .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
795: @*/
796: PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
797: {
798: DM_Forest *forest = (DM_Forest *) dm->data;
803: *initRefinement = forest->initRefinement;
804: return(0);
805: }
809: /*@
810: DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
811: DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest
812: (see DMForestGetAdaptivityForest()), this limits the amount of refinement.
814: Logically collective on dm
816: Input Parameters:
817: + dm - the forest
818: - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
820: Level: intermediate
822: .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM()
823: @*/
824: PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
825: {
826: DM_Forest *forest = (DM_Forest *) dm->data;
830: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup");
831: forest->maxRefinement = maxRefinement;
832: return(0);
833: }
837: /*@
838: DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
839: DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest (see
840: DMForestGetAdaptivityForest()), this limits the amount of refinement.
842: Not collective
844: Input Parameter:
845: . dm - the forest
847: Output Parameter:
848: . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
850: Level: intermediate
852: .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
853: @*/
854: PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
855: {
856: DM_Forest *forest = (DM_Forest *) dm->data;
861: *maxRefinement = forest->maxRefinement;
862: return(0);
863: }
867: /*@C
868: DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
869: Subtypes of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
870: for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
871: specify refinement/coarsening.
873: Logically collective on dm
875: Input Parameters:
876: + dm - the forest
877: - adaptStrategy - default DMFORESTADAPTALL
879: Level: advanced
881: .seealso: DMForestGetAdaptivityStrategy()
882: @*/
883: PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
884: {
885: DM_Forest *forest = (DM_Forest *) dm->data;
890: PetscFree(forest->adaptStrategy);
891: PetscStrallocpy((const char *) adaptStrategy,(char **)&forest->adaptStrategy);
892: return(0);
893: }
897: /*@C
898: DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes. Subtypes
899: of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all
900: processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
901: one process needs to specify refinement/coarsening.
903: Not collective
905: Input Parameter:
906: . dm - the forest
908: Output Parameter:
909: . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
911: Level: advanced
913: .seealso: DMForestSetAdaptivityStrategy()
914: @*/
915: PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
916: {
917: DM_Forest *forest = (DM_Forest *) dm->data;
922: *adaptStrategy = forest->adaptStrategy;
923: return(0);
924: }
928: /*@
929: DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
930: relating the cells of the pre-adaptation forest to the post-adaptiation forest. After DMSetUp() is called, these transfer PetscSFs can be accessed with DMForestGetAdaptivitySF().
932: Logically collective on dm
934: Input Parameters:
935: + dm - the post-adaptation forest
936: - computeSF - default PETSC_TRUE
938: Level: advanced
940: .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
941: @*/
942: PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
943: {
944: DM_Forest *forest;
948: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called");
949: forest = (DM_Forest *) dm->data;
950: forest->computeAdaptSF = computeSF;
951: return(0);
952: }
956: /*@
957: DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
958: pre-adaptation forest to the post-adaptiation forest. After DMSetUp() is called, these transfer PetscSFs can be
959: accessed with DMForestGetAdaptivitySF().
961: Not collective
963: Input Parameter:
964: . dm - the post-adaptation forest
966: Output Parameter:
967: . computeSF - default PETSC_TRUE
969: Level: advanced
971: .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
972: @*/
973: PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
974: {
975: DM_Forest *forest;
979: forest = (DM_Forest *) dm->data;
980: *computeSF = forest->computeAdaptSF;
981: return(0);
982: }
986: /*@
987: DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
988: Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
989: some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa. Therefore there are two
990: PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
991: pre-adaptation fine cells to post-adaptation coarse cells.
993: Not collective
995: Input Parameter:
996: dm - the post-adaptation forest
998: Output Parameter:
999: preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
1000: coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1002: Level: advanced
1004: .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF()
1005: @*/
1006: PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1007: {
1008: DM_Forest *forest;
1013: DMSetUp(dm);
1014: forest = (DM_Forest *) dm->data;
1015: if (preCoarseToFine) {
1016: *preCoarseToFine = forest->preCoarseToFine;
1017: }
1018: if (coarseToPreFine) {
1019: *coarseToPreFine = forest->coarseToPreFine;
1020: }
1021: return(0);
1022: }
1026: /*@
1027: DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
1028: indicate that the diameter of neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may
1029: only support one particular choice of grading factor.
1031: Logically collective on dm
1033: Input Parameters:
1034: + dm - the forest
1035: - grade - the grading factor
1037: Level: advanced
1039: .seealso: DMForestGetGradeFactor()
1040: @*/
1041: PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1042: {
1043: DM_Forest *forest = (DM_Forest *) dm->data;
1047: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1048: forest->gradeFactor = grade;
1049: return(0);
1050: }
1054: /*@
1055: DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
1056: neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may only support one particular
1057: choice of grading factor.
1059: Not collective
1061: Input Parameter:
1062: . dm - the forest
1064: Output Parameter:
1065: . grade - the grading factor
1067: Level: advanced
1069: .seealso: DMForestSetGradeFactor()
1070: @*/
1071: PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1072: {
1073: DM_Forest *forest = (DM_Forest *) dm->data;
1078: *grade = forest->gradeFactor;
1079: return(0);
1080: }
1084: /*@
1085: DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
1086: the cell weight (see DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be
1087: (cellWeight) * (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on
1088: its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
1089: computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
1091: Logically collective on dm
1093: Input Parameters:
1094: + dm - the forest
1095: - weightsFactors - default 1.
1097: Level: advanced
1099: .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
1100: @*/
1101: PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1102: {
1103: DM_Forest *forest = (DM_Forest *) dm->data;
1107: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1108: forest->weightsFactor = weightsFactor;
1109: return(0);
1110: }
1114: /*@
1115: DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
1116: DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be (cellWeight) *
1117: (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on its level; a
1118: factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
1119: associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
1121: Not collective
1123: Input Parameter:
1124: . dm - the forest
1126: Output Parameter:
1127: . weightsFactors - default 1.
1129: Level: advanced
1131: .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights()
1132: @*/
1133: PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1134: {
1135: DM_Forest *forest = (DM_Forest *) dm->data;
1140: *weightsFactor = forest->weightsFactor;
1141: return(0);
1142: }
1146: /*@
1147: DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
1149: Not collective
1151: Input Parameter:
1152: . dm - the forest
1154: Output Parameters:
1155: + cStart - the first cell on this process
1156: - cEnd - one after the final cell on this process
1158: Level: intermediate
1160: .seealso: DMForestGetCellSF()
1161: @*/
1162: PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1163: {
1164: DM_Forest *forest = (DM_Forest *) dm->data;
1171: if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1172: forest->createcellchart(dm,&forest->cStart,&forest->cEnd);
1173: }
1174: *cStart = forest->cStart;
1175: *cEnd = forest->cEnd;
1176: return(0);
1177: }
1181: /*@
1182: DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
1184: Not collective
1186: Input Parameter:
1187: . dm - the forest
1189: Output Parameter:
1190: . cellSF - the PetscSF
1192: Level: intermediate
1194: .seealso: DMForestGetCellChart()
1195: @*/
1196: PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1197: {
1198: DM_Forest *forest = (DM_Forest *) dm->data;
1204: if ((!forest->cellSF) && forest->createcellsf) {
1205: forest->createcellsf(dm,&forest->cellSF);
1206: }
1207: *cellSF = forest->cellSF;
1208: return(0);
1209: }
1213: /*@C
1214: DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
1215: DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination). The
1216: interpretation of the label values is up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and
1217: DM_FOREST_COARSEN have been reserved as choices that should be accepted by all subtypes.
1219: Logically collective on dm
1221: Input Parameters:
1222: - dm - the forest
1223: + adaptLabel - the name of the label in the pre-adaptation forest
1225: Level: intermediate
1227: .seealso DMForestGetAdaptivityLabel()
1228: @*/
1229: PetscErrorCode DMForestSetAdaptivityLabel(DM dm, const char * adaptLabel)
1230: {
1231: DM_Forest *forest = (DM_Forest *) dm->data;
1236: PetscFree(forest->adaptLabel);
1237: PetscStrallocpy(adaptLabel,&forest->adaptLabel);
1238: return(0);
1239: }
1243: /*@C
1244: DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
1245: holds the adaptation flags (refinement, coarsening, or some combination). The interpretation of the label values is
1246: up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and DM_FOREST_COARSEN have been reserved as
1247: choices that should be accepted by all subtypes.
1249: Not collective
1251: Input Parameter:
1252: . dm - the forest
1254: Output Parameter:
1255: . adaptLabel - the name of the label in the pre-adaptation forest
1257: Level: intermediate
1259: .seealso DMForestSetAdaptivityLabel()
1260: @*/
1261: PetscErrorCode DMForestGetAdaptivityLabel(DM dm, const char ** adaptLabel)
1262: {
1263: DM_Forest *forest = (DM_Forest *) dm->data;
1267: *adaptLabel = forest->adaptLabel;
1268: return(0);
1269: }
1273: /*@
1274: DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1275: process: weights are used to determine parallel partitioning. Partitions will be created so that each process's
1276: ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1277: of 1.
1279: Logically collective on dm
1281: Input Parameters:
1282: + dm - the forest
1283: . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1284: - copyMode - how weights should reference weights
1286: Level: advanced
1288: .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
1289: @*/
1290: PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1291: {
1292: DM_Forest *forest = (DM_Forest *) dm->data;
1293: PetscInt cStart, cEnd;
1298: DMForestGetCellChart(dm,&cStart,&cEnd);
1299: if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1300: if (copyMode == PETSC_COPY_VALUES) {
1301: if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1302: PetscMalloc1(cEnd-cStart,&forest->cellWeights);
1303: }
1304: PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));
1305: forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1306: return(0);
1307: }
1308: if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1309: PetscFree(forest->cellWeights);
1310: }
1311: forest->cellWeights = weights;
1312: forest->cellWeightsCopyMode = copyMode;
1313: return(0);
1314: }
1318: /*@
1319: DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1320: process: weights are used to determine parallel partitioning. Partitions will be created so that each process's
1321: ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1322: of 1.
1324: Not collective
1326: Input Parameter:
1327: . dm - the forest
1329: Output Parameter:
1330: . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1332: Level: advanced
1334: .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
1335: @*/
1336: PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1337: {
1338: DM_Forest *forest = (DM_Forest *) dm->data;
1343: *weights = forest->cellWeights;
1344: return(0);
1345: }
1349: /*@
1350: DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
1351: a pre-adaptation forest (see DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each
1352: process's cells to the process's capacity will be roughly equal for all processes. A capacity of 0 indicates that
1353: the current process should not have any cells after repartitioning.
1355: Logically Collective on dm
1357: Input parameters:
1358: + dm - the forest
1359: - capacity - this process's capacity
1361: Level: advanced
1363: .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1364: @*/
1365: PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1366: {
1367: DM_Forest *forest = (DM_Forest *) dm->data;
1371: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1372: if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1373: forest->weightCapacity = capacity;
1374: return(0);
1375: }
1379: /*@
1380: DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
1381: DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each process's cells to the
1382: process's capacity will be roughly equal for all processes. A capacity of 0 indicates that the current process
1383: should not have any cells after repartitioning.
1385: Not collective
1387: Input parameter:
1388: . dm - the forest
1390: Output parameter:
1391: . capacity - this process's capacity
1393: Level: advanced
1395: .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1396: @*/
1397: PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1398: {
1399: DM_Forest *forest = (DM_Forest *) dm->data;
1404: *capacity = forest->weightCapacity;
1405: return(0);
1406: }
1410: PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1411: {
1412: DM_Forest *forest = (DM_Forest *) dm->data;
1413: PetscBool flg, flg1, flg2, flg3, flg4;
1414: DMForestTopology oldTopo;
1415: char stringBuffer[256];
1416: PetscViewer viewer;
1417: PetscViewerFormat format;
1418: PetscInt adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1419: PetscReal weightsFactor;
1420: DMForestAdaptivityStrategy adaptStrategy;
1421: PetscErrorCode ierr;
1425: forest->setfromoptionscalled = PETSC_TRUE;
1426: DMForestGetTopology(dm, &oldTopo);
1427: PetscOptionsHead(PetscOptionsObject,"DMForest Options");
1428: PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);
1429: PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);
1430: PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);
1431: PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);
1432: if ((PetscInt) flg1 + (PetscInt) flg2 + (PetscInt) flg3 + (PetscInt) flg4 > 1) {
1433: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}");
1434: }
1435: if (flg1) {
1436: DMForestSetTopology(dm,(DMForestTopology)stringBuffer);
1437: DMForestSetBaseDM(dm,NULL);
1438: DMForestSetAdaptivityForest(dm,NULL);
1439: }
1440: if (flg2) {
1441: DM base;
1443: DMCreate(PetscObjectComm((PetscObject)dm),&base);
1444: PetscViewerPushFormat(viewer,format);
1445: DMLoad(base,viewer);
1446: PetscViewerDestroy(&viewer);
1447: DMForestSetBaseDM(dm,base);
1448: DMDestroy(&base);
1449: DMForestSetTopology(dm,NULL);
1450: DMForestSetAdaptivityForest(dm,NULL);
1451: }
1452: if (flg3) {
1453: DM coarse;
1455: DMCreate(PetscObjectComm((PetscObject)dm),&coarse);
1456: PetscViewerPushFormat(viewer,format);
1457: DMLoad(coarse,viewer);
1458: PetscViewerDestroy(&viewer);
1459: DMForestSetAdaptivityForest(dm,coarse);
1460: DMDestroy(&coarse);
1461: DMForestSetTopology(dm,NULL);
1462: DMForestSetBaseDM(dm,NULL);
1463: }
1464: if (flg4) {
1465: DM fine;
1467: DMCreate(PetscObjectComm((PetscObject)dm),&fine);
1468: PetscViewerPushFormat(viewer,format);
1469: DMLoad(fine,viewer);
1470: PetscViewerDestroy(&viewer);
1471: DMForestSetAdaptivityForest(dm,fine);
1472: DMDestroy(&fine);
1473: DMForestSetTopology(dm,NULL);
1474: DMForestSetBaseDM(dm,NULL);
1475: }
1476: DMForestGetAdjacencyDimension(dm,&adjDim);
1477: PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);
1478: if (flg) {
1479: DMForestSetAdjacencyDimension(dm,adjDim);
1480: }
1481: else {
1482: DMForestGetAdjacencyCodimension(dm,&adjCodim);
1483: PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);
1484: if (flg) {
1485: DMForestSetAdjacencyCodimension(dm,adjCodim);
1486: }
1487: }
1488: DMForestGetPartitionOverlap(dm,&overlap);
1489: PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);
1490: if (flg) {
1491: DMForestSetPartitionOverlap(dm,overlap);
1492: }
1493: #if 0
1494: PetscOptionsInt("-dm_refine","equivalent to -dm_forest_set_minimum_refinement and -dm_forest_set_initial_refinement with the same value",NULL,minRefinement,&minRefinement,&flg);
1495: if (flg) {
1496: DMForestSetMinimumRefinement(dm,minRefinement);
1497: DMForestSetInitialRefinement(dm,minRefinement);
1498: }
1499: PetscOptionsInt("-dm_refine_hierarchy","equivalent to -dm_forest_set_minimum_refinement 0 and -dm_forest_set_initial_refinement",NULL,initRefinement,&initRefinement,&flg);
1500: if (flg) {
1501: DMForestSetMinimumRefinement(dm,0);
1502: DMForestSetInitialRefinement(dm,initRefinement);
1503: }
1504: #endif
1505: DMForestGetMinimumRefinement(dm,&minRefinement);
1506: PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);
1507: if (flg) {
1508: DMForestSetMinimumRefinement(dm,minRefinement);
1509: }
1510: DMForestGetInitialRefinement(dm,&initRefinement);
1511: PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);
1512: if (flg) {
1513: DMForestSetInitialRefinement(dm,initRefinement);
1514: }
1515: DMForestGetMaximumRefinement(dm,&maxRefinement);
1516: PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);
1517: if (flg) {
1518: DMForestSetMaximumRefinement(dm,maxRefinement);
1519: }
1520: DMForestGetAdaptivityStrategy(dm,&adaptStrategy);
1521: PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);
1522: if (flg) {
1523: DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);
1524: }
1525: DMForestGetGradeFactor(dm,&grade);
1526: PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);
1527: if (flg) {
1528: DMForestSetGradeFactor(dm,grade);
1529: }
1530: DMForestGetCellWeightFactor(dm,&weightsFactor);
1531: PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);
1532: if (flg) {
1533: DMForestSetCellWeightFactor(dm,weightsFactor);
1534: }
1535: PetscOptionsTail();
1536: return(0);
1537: }
1541: PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1542: {
1546: if (subdm) {DMClone(dm, subdm);}
1547: DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);
1548: return(0);
1549: }
1553: PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
1554: {
1555: DMLabel refine;
1556: DM fineDM;
1560: DMGetFineDM(dm,&fineDM);
1561: if (fineDM) {
1562: PetscObjectReference((PetscObject)fineDM);
1563: *dmRefined = fineDM;
1564: return(0);
1565: }
1566: DMForestTemplate(dm,comm,dmRefined);
1567: DMGetLabel(dm,"refine",&refine);
1568: if (!refine) {
1569: DMCreateLabel(dm,"refine");
1570: DMGetLabel(dm,"refine",&refine);
1571: DMLabelSetDefaultValue(refine,DM_FOREST_REFINE);
1572: }
1573: DMForestSetAdaptivityLabel(*dmRefined,"refine");
1574: return(0);
1575: }
1579: PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
1580: {
1581: DMLabel coarsen;
1582: DM coarseDM;
1586: {
1587: PetscMPIInt mpiComparison;
1588: MPI_Comm dmcomm = PetscObjectComm((PetscObject)dm);
1590: MPI_Comm_compare(comm, dmcomm, &mpiComparison);
1591: if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) {
1592: SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
1593: }
1594: }
1595: DMGetCoarseDM(dm,&coarseDM);
1596: if (coarseDM) {
1597: PetscObjectReference((PetscObject)coarseDM);
1598: *dmCoarsened = coarseDM;
1599: return(0);
1600: }
1601: DMForestTemplate(dm,comm,dmCoarsened);
1602: DMForestSetAdaptivityPurpose(coarseDM,DM_FOREST_COARSEN);
1603: DMGetLabel(dm,"coarsen",&coarsen);
1604: if (!coarsen) {
1605: DMCreateLabel(dm,"coarsen");
1606: DMGetLabel(dm,"coarsen",&coarsen);
1607: DMLabelSetDefaultValue(coarsen,DM_FOREST_COARSEN);
1608: }
1609: DMForestSetAdaptivityLabel(*dmCoarsened,"coarsen");
1610: return(0);
1611: }
1615: static PetscErrorCode DMInitialize_Forest(DM dm)
1616: {
1620: PetscMemzero(dm->ops,sizeof(*(dm->ops)));
1622: dm->ops->clone = DMClone_Forest;
1623: dm->ops->setfromoptions = DMSetFromOptions_Forest;
1624: dm->ops->destroy = DMDestroy_Forest;
1625: dm->ops->createsubdm = DMCreateSubDM_Forest;
1626: dm->ops->refine = DMRefine_Forest;
1627: dm->ops->coarsen = DMCoarsen_Forest;
1628: return(0);
1629: }
1631: /*MC
1632: DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh. Forests usually have a base DM (see DMForestGetBaseDM()), from which it is refined. The refinement and partitioning of forests is considered immutable after DMSetUp() is called. To adapt a mesh, one should call DMForestTemplate() to create a new mesh that will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new mesh.
1634: To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the previous mesh whose values indicate which cells should be refined (DM_FOREST_REFINE) or coarsened (DM_FOREST_COARSEN) and how (subtypes are free to allow additional values for things like anisotropic refinement). The name of the label should be given to the *new* mesh with DMForestSetAdaptivityLabel().
1636: Level: advanced
1638: .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
1639: M*/
1643: PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1644: {
1645: DM_Forest *forest;
1650: PetscNewLog(dm,&forest);
1651: dm->dim = 0;
1652: dm->data = forest;
1653: forest->refct = 1;
1654: forest->data = NULL;
1655: forest->setfromoptionscalled = PETSC_FALSE;
1656: forest->topology = NULL;
1657: forest->base = NULL;
1658: forest->adjDim = PETSC_DEFAULT;
1659: forest->overlap = PETSC_DEFAULT;
1660: forest->minRefinement = PETSC_DEFAULT;
1661: forest->maxRefinement = PETSC_DEFAULT;
1662: forest->initRefinement = PETSC_DEFAULT;
1663: forest->cStart = PETSC_DETERMINE;
1664: forest->cEnd = PETSC_DETERMINE;
1665: forest->cellSF = 0;
1666: forest->adaptLabel = NULL;
1667: forest->gradeFactor = 2;
1668: forest->cellWeights = NULL;
1669: forest->cellWeightsCopyMode = PETSC_USE_POINTER;
1670: forest->weightsFactor = 1.;
1671: forest->weightCapacity = 1.;
1672: DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);
1673: DMInitialize_Forest(dm);
1674: return(0);
1675: }