xref: /petsc/src/dm/dt/dualspace/impls/refined/dualspacerefined.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
11ac17e89SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
21ac17e89SToby Isaac #include <petscdmplex.h>
31ac17e89SToby Isaac 
41ac17e89SToby Isaac typedef struct {
51ac17e89SToby Isaac   PetscInt dummy;
61ac17e89SToby Isaac } PetscDualSpace_Refined;
71ac17e89SToby Isaac 
81ac17e89SToby Isaac /*@
91ac17e89SToby Isaac   PetscDualSpaceRefinedSetCellSpaces - Set the dual spaces for the closures of each of the cells
10dce8aebaSBarry Smith   in the multicell `DM` of a `PetscDualSpace`
111ac17e89SToby Isaac 
1220f4b53cSBarry Smith   Collective
131ac17e89SToby Isaac 
144165533cSJose E. Roman   Input Parameters:
15dce8aebaSBarry Smith + sp         - a `PetscDualSpace`
16dce8aebaSBarry Smith - cellSpaces - one `PetscDualSpace` for each of the cells.  The reference count of each cell space will be incremented,
171ac17e89SToby Isaac                 so the user is still responsible for these spaces afterwards
181ac17e89SToby Isaac 
191ac17e89SToby Isaac   Level: intermediate
201ac17e89SToby Isaac 
21b24fb147SBarry Smith .seealso: `PETSCDUALSPACEREFINED`, `PetscDualSpace`, `PetscFERefine()`
221ac17e89SToby Isaac @*/
PetscDualSpaceRefinedSetCellSpaces(PetscDualSpace sp,const PetscDualSpace cellSpaces[])23d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDualSpaceRefinedSetCellSpaces(PetscDualSpace sp, const PetscDualSpace cellSpaces[])
24d71ae5a4SJacob Faibussowitsch {
251ac17e89SToby Isaac   PetscFunctionBegin;
261ac17e89SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
274f572ea9SToby Isaac   PetscAssertPointer(cellSpaces, 2);
2828b400f6SJacob Faibussowitsch   PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change cell spaces after setup is called");
29cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceRefinedSetCellSpaces_C", (PetscDualSpace, const PetscDualSpace[]), (sp, cellSpaces));
303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
311ac17e89SToby Isaac }
321ac17e89SToby Isaac 
PetscDualSpaceRefinedSetCellSpaces_Refined(PetscDualSpace sp,const PetscDualSpace cellSpaces[])33d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceRefinedSetCellSpaces_Refined(PetscDualSpace sp, const PetscDualSpace cellSpaces[])
34d71ae5a4SJacob Faibussowitsch {
351ac17e89SToby Isaac   DM       dm;
361ac17e89SToby Isaac   PetscInt pStart, pEnd;
371ac17e89SToby Isaac   PetscInt cStart, cEnd, c;
381ac17e89SToby Isaac 
391ac17e89SToby Isaac   PetscFunctionBegin;
401ac17e89SToby Isaac   dm = sp->dm;
4128b400f6SJacob Faibussowitsch   PetscCheck(dm, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "PetscDualSpace must have a DM (PetscDualSpaceSetDM()) before calling PetscDualSpaceRefinedSetCellSpaces");
429566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
43f4f49eeaSPierre Jolivet   if (!sp->pointSpaces) PetscCall(PetscCalloc1(pEnd - pStart, &sp->pointSpaces));
449566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
451ac17e89SToby Isaac   for (c = 0; c < cEnd - cStart; c++) {
469566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)cellSpaces[c]));
47f4f49eeaSPierre Jolivet     PetscCall(PetscDualSpaceDestroy(&sp->pointSpaces[c + cStart - pStart]));
481ac17e89SToby Isaac     sp->pointSpaces[c + cStart - pStart] = cellSpaces[c];
491ac17e89SToby Isaac   }
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
511ac17e89SToby Isaac }
521ac17e89SToby Isaac 
PetscDualSpaceDestroy_Refined(PetscDualSpace sp)53d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceDestroy_Refined(PetscDualSpace sp)
54d71ae5a4SJacob Faibussowitsch {
551ac17e89SToby Isaac   PetscDualSpace_Refined *ref = (PetscDualSpace_Refined *)sp->data;
561ac17e89SToby Isaac 
571ac17e89SToby Isaac   PetscFunctionBegin;
589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceRefinedSetCellSpaces_C", NULL));
599566063dSJacob Faibussowitsch   PetscCall(PetscFree(ref));
603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
611ac17e89SToby Isaac }
621ac17e89SToby Isaac 
PetscDualSpaceSetUp_Refined(PetscDualSpace sp)63d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceSetUp_Refined(PetscDualSpace sp)
64d71ae5a4SJacob Faibussowitsch {
651ac17e89SToby Isaac   PetscInt     pStart, pEnd, depth;
661ac17e89SToby Isaac   PetscInt     cStart, cEnd, c, spdim;
671ac17e89SToby Isaac   PetscInt     h;
681ac17e89SToby Isaac   DM           dm;
691ac17e89SToby Isaac   PetscSection section;
701ac17e89SToby Isaac 
711ac17e89SToby Isaac   PetscFunctionBegin;
729566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(sp, &dm));
739566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
749566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
759566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
761ac17e89SToby Isaac   for (c = cStart; c < cEnd; c++) {
771ac17e89SToby Isaac     if (sp->pointSpaces[c - pStart]) {
781ac17e89SToby Isaac       PetscInt ccStart, ccEnd;
7908401ef6SPierre Jolivet       PetscCheck(sp->pointSpaces[c - pStart]->k == sp->k, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "All cell spaces must have the same form degree as the refined dual space");
8008401ef6SPierre Jolivet       PetscCheck(sp->pointSpaces[c - pStart]->Nc == sp->Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "All cell spaces must have the same number of components as the refined dual space");
819566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(sp->pointSpaces[c - pStart]->dm, 0, &ccStart, &ccEnd));
821dca8a05SBarry Smith       PetscCheck(ccEnd - ccStart == 1, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "All cell spaces must have a single cell themselves");
831ac17e89SToby Isaac     }
841ac17e89SToby Isaac   }
851ac17e89SToby Isaac   for (c = cStart; c < cEnd; c++) {
861ac17e89SToby Isaac     if (sp->pointSpaces[c - pStart]) {
871ac17e89SToby Isaac       PetscBool cUniform;
881ac17e89SToby Isaac 
899566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetUniform(sp->pointSpaces[c - pStart], &cUniform));
901ac17e89SToby Isaac       if (!cUniform) break;
911ac17e89SToby Isaac     }
921ac17e89SToby Isaac     if ((c > cStart) && sp->pointSpaces[c - pStart] != sp->pointSpaces[c - 1 - pStart]) break;
931ac17e89SToby Isaac   }
941ac17e89SToby Isaac   if (c < cEnd) sp->uniform = PETSC_FALSE;
951ac17e89SToby Isaac   for (h = 0; h < depth; h++) {
961ac17e89SToby Isaac     PetscInt hStart, hEnd;
971ac17e89SToby Isaac 
989566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd));
991ac17e89SToby Isaac     for (c = hStart; c < hEnd; c++) {
1001ac17e89SToby Isaac       PetscInt        coneSize, e;
1011ac17e89SToby Isaac       PetscDualSpace  cspace = sp->pointSpaces[c - pStart];
1021ac17e89SToby Isaac       const PetscInt *cone;
1031ac17e89SToby Isaac       const PetscInt *refCone;
1041ac17e89SToby Isaac 
1051ac17e89SToby Isaac       if (!cspace) continue;
1069566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, c, &coneSize));
1079566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, c, &cone));
1089566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(cspace->dm, 0, &refCone));
1091ac17e89SToby Isaac       for (e = 0; e < coneSize; e++) {
1101ac17e89SToby Isaac         PetscInt       point    = cone[e];
1111ac17e89SToby Isaac         PetscInt       refpoint = refCone[e];
1121ac17e89SToby Isaac         PetscDualSpace espace;
1131ac17e89SToby Isaac 
1149566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetPointSubspace(cspace, refpoint, &espace));
1151ac17e89SToby Isaac         if (sp->pointSpaces[point - pStart] == NULL) {
1169566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)espace));
1171ac17e89SToby Isaac           sp->pointSpaces[point - pStart] = espace;
1181ac17e89SToby Isaac         }
1191ac17e89SToby Isaac       }
1201ac17e89SToby Isaac     }
1211ac17e89SToby Isaac   }
1229566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetSection(sp, &section));
1239566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(sp, &spdim));
124f4f49eeaSPierre Jolivet   PetscCall(PetscMalloc1(spdim, &sp->functional));
1259566063dSJacob Faibussowitsch   PetscCall(PetscDualSpacePushForwardSubspaces_Internal(sp, pStart, pEnd));
1263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1271ac17e89SToby Isaac }
1281ac17e89SToby Isaac 
PetscDualSpaceRefinedView_Ascii(PetscDualSpace sp,PetscViewer viewer)129d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceRefinedView_Ascii(PetscDualSpace sp, PetscViewer viewer)
130d71ae5a4SJacob Faibussowitsch {
1311ac17e89SToby Isaac   PetscFunctionBegin;
1321ac17e89SToby Isaac   if (sp->dm && sp->pointSpaces) {
1331ac17e89SToby Isaac     PetscInt pStart, pEnd;
1341ac17e89SToby Isaac     PetscInt cStart, cEnd, c;
1351ac17e89SToby Isaac 
1369566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(sp->dm, &pStart, &pEnd));
1379566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(sp->dm, 0, &cStart, &cEnd));
1389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Refined dual space:\n"));
1399566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1401ac17e89SToby Isaac     for (c = cStart; c < cEnd; c++) {
1411ac17e89SToby Isaac       if (!sp->pointSpaces[c - pStart]) {
14263a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Cell space %" PetscInt_FMT " not set yet\n", c));
1431ac17e89SToby Isaac       } else {
14463a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Cell space %" PetscInt_FMT ":ot set yet\n", c));
1459566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
1469566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceView(sp->pointSpaces[c - pStart], viewer));
1479566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
1481ac17e89SToby Isaac       }
1491ac17e89SToby Isaac     }
1509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1511baa6e33SBarry Smith   } else PetscCall(PetscViewerASCIIPrintf(viewer, "Refined dual space: (cell spaces not set yet)\n"));
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1531ac17e89SToby Isaac }
1541ac17e89SToby Isaac 
PetscDualSpaceView_Refined(PetscDualSpace sp,PetscViewer viewer)155d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceView_Refined(PetscDualSpace sp, PetscViewer viewer)
156d71ae5a4SJacob Faibussowitsch {
157*9f196a02SMartin Diehl   PetscBool isascii;
1581ac17e89SToby Isaac 
1591ac17e89SToby Isaac   PetscFunctionBegin;
1601ac17e89SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1611ac17e89SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
162*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
163*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscDualSpaceRefinedView_Ascii(sp, viewer));
1643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1651ac17e89SToby Isaac }
1661ac17e89SToby Isaac 
PetscDualSpaceInitialize_Refined(PetscDualSpace sp)167d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceInitialize_Refined(PetscDualSpace sp)
168d71ae5a4SJacob Faibussowitsch {
1691ac17e89SToby Isaac   PetscFunctionBegin;
1701ac17e89SToby Isaac   sp->ops->destroy              = PetscDualSpaceDestroy_Refined;
1711ac17e89SToby Isaac   sp->ops->view                 = PetscDualSpaceView_Refined;
1721ac17e89SToby Isaac   sp->ops->setfromoptions       = NULL;
1731ac17e89SToby Isaac   sp->ops->duplicate            = NULL;
1741ac17e89SToby Isaac   sp->ops->setup                = PetscDualSpaceSetUp_Refined;
1751ac17e89SToby Isaac   sp->ops->createheightsubspace = NULL;
1761ac17e89SToby Isaac   sp->ops->createpointsubspace  = NULL;
1771ac17e89SToby Isaac   sp->ops->getsymmetries        = NULL;
1781ac17e89SToby Isaac   sp->ops->apply                = PetscDualSpaceApplyDefault;
1791ac17e89SToby Isaac   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
1801ac17e89SToby Isaac   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
1811ac17e89SToby Isaac   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
1821ac17e89SToby Isaac   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1841ac17e89SToby Isaac }
1851ac17e89SToby Isaac 
1861ac17e89SToby Isaac /*MC
187dce8aebaSBarry Smith   PETSCDUALSPACEREFINED = "refined" - A `PetscDualSpaceType` that defines the joint dual space of a group of cells, usually refined from one larger cell
1881ac17e89SToby Isaac 
1891ac17e89SToby Isaac   Level: intermediate
1901ac17e89SToby Isaac 
191b24fb147SBarry Smith .seealso: `PetscDualSpace`, `PetscDualSpaceType`, `PetscDualSpaceRefinedSetCellSpaces`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`
1921ac17e89SToby Isaac M*/
PetscDualSpaceCreate_Refined(PetscDualSpace sp)193d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Refined(PetscDualSpace sp)
194d71ae5a4SJacob Faibussowitsch {
1951ac17e89SToby Isaac   PetscDualSpace_Refined *ref;
1961ac17e89SToby Isaac 
1971ac17e89SToby Isaac   PetscFunctionBegin;
1981ac17e89SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1994dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ref));
2001ac17e89SToby Isaac   sp->data = ref;
2011ac17e89SToby Isaac 
2029566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceInitialize_Refined(sp));
2039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceRefinedSetCellSpaces_C", PetscDualSpaceRefinedSetCellSpaces_Refined));
2043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2051ac17e89SToby Isaac }
206