xref: /petsc/src/dm/dt/dualspace/impls/refined/dualspacerefined.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
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
101ac17e89SToby Isaac    in the multicell DM of a PetscDualSpace
111ac17e89SToby Isaac 
121ac17e89SToby Isaac    Collective on PetscDualSpace
131ac17e89SToby Isaac 
144165533cSJose E. Roman    Input Parameters:
151ac17e89SToby Isaac +  sp - a PetscDualSpace
161ac17e89SToby Isaac -  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 
211ac17e89SToby Isaac .seealso: PetscFERefine()
221ac17e89SToby Isaac @*/
231ac17e89SToby Isaac PetscErrorCode PetscDualSpaceRefinedSetCellSpaces(PetscDualSpace sp, const PetscDualSpace cellSpaces[])
241ac17e89SToby Isaac {
251ac17e89SToby Isaac   PetscFunctionBegin;
261ac17e89SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
271ac17e89SToby Isaac   PetscValidPointer(cellSpaces,2);
2828b400f6SJacob Faibussowitsch   PetscCheck(!sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change cell spaces after setup is called");
29*cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceRefinedSetCellSpaces_C", (PetscDualSpace,const PetscDualSpace []),(sp,cellSpaces));
301ac17e89SToby Isaac   PetscFunctionReturn(0);
311ac17e89SToby Isaac }
321ac17e89SToby Isaac 
331ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceRefinedSetCellSpaces_Refined(PetscDualSpace sp, const PetscDualSpace cellSpaces[])
341ac17e89SToby Isaac {
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));
431ac17e89SToby Isaac   if (!sp->pointSpaces) {
441ac17e89SToby Isaac 
459566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(pEnd-pStart,&(sp->pointSpaces)));
461ac17e89SToby Isaac   }
479566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
481ac17e89SToby Isaac   for (c = 0; c < cEnd - cStart; c++) {
499566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)cellSpaces[c]));
509566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceDestroy(&(sp->pointSpaces[c + cStart - pStart])));
511ac17e89SToby Isaac     sp->pointSpaces[c+cStart-pStart] = cellSpaces[c];
521ac17e89SToby Isaac   }
531ac17e89SToby Isaac   PetscFunctionReturn(0);
541ac17e89SToby Isaac }
551ac17e89SToby Isaac 
561ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceDestroy_Refined(PetscDualSpace sp)
571ac17e89SToby Isaac {
581ac17e89SToby Isaac   PetscDualSpace_Refined *ref = (PetscDualSpace_Refined *) sp->data;
591ac17e89SToby Isaac 
601ac17e89SToby Isaac   PetscFunctionBegin;
619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceRefinedSetCellSpaces_C", NULL));
629566063dSJacob Faibussowitsch   PetscCall(PetscFree(ref));
631ac17e89SToby Isaac   PetscFunctionReturn(0);
641ac17e89SToby Isaac }
651ac17e89SToby Isaac 
661ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceSetUp_Refined(PetscDualSpace sp)
671ac17e89SToby Isaac {
681ac17e89SToby Isaac   PetscInt pStart, pEnd, depth;
691ac17e89SToby Isaac   PetscInt cStart, cEnd, c, spdim;
701ac17e89SToby Isaac   PetscInt h;
711ac17e89SToby Isaac   DM dm;
721ac17e89SToby Isaac   PetscSection   section;
731ac17e89SToby Isaac 
741ac17e89SToby Isaac   PetscFunctionBegin;
759566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(sp, &dm));
769566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
779566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
789566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
791ac17e89SToby Isaac   for (c = cStart; c < cEnd; c++) {
801ac17e89SToby Isaac     if (sp->pointSpaces[c-pStart]) {
811ac17e89SToby Isaac       PetscInt ccStart, ccEnd;
822c71b3e2SJacob Faibussowitsch       PetscCheckFalse(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");
832c71b3e2SJacob Faibussowitsch       PetscCheckFalse(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");
849566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(sp->pointSpaces[c-pStart]->dm, 0, &ccStart, &ccEnd));
852c71b3e2SJacob Faibussowitsch       PetscCheckFalse(ccEnd - ccStart != 1,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "All cell spaces must have a single cell themselves");
861ac17e89SToby Isaac     }
871ac17e89SToby Isaac   }
881ac17e89SToby Isaac   for (c = cStart; c < cEnd; c++) {
891ac17e89SToby Isaac     if (sp->pointSpaces[c-pStart]) {
901ac17e89SToby Isaac       PetscBool cUniform;
911ac17e89SToby Isaac 
929566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetUniform(sp->pointSpaces[c-pStart],&cUniform));
931ac17e89SToby Isaac       if (!cUniform) break;
941ac17e89SToby Isaac     }
951ac17e89SToby Isaac     if ((c > cStart) && sp->pointSpaces[c-pStart] != sp->pointSpaces[c-1-pStart]) break;
961ac17e89SToby Isaac   }
971ac17e89SToby Isaac   if (c < cEnd) sp->uniform = PETSC_FALSE;
981ac17e89SToby Isaac   for (h = 0; h < depth; h++) {
991ac17e89SToby Isaac     PetscInt hStart, hEnd;
1001ac17e89SToby Isaac 
1019566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd));
1021ac17e89SToby Isaac     for (c = hStart; c < hEnd; c++) {
1031ac17e89SToby Isaac       PetscInt coneSize, e;
1041ac17e89SToby Isaac       PetscDualSpace cspace = sp->pointSpaces[c-pStart];
1051ac17e89SToby Isaac       const PetscInt *cone;
1061ac17e89SToby Isaac       const PetscInt *refCone;
1071ac17e89SToby Isaac 
1081ac17e89SToby Isaac       if (!cspace) continue;
1099566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, c, &coneSize));
1109566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, c, &cone));
1119566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(cspace->dm, 0, &refCone));
1121ac17e89SToby Isaac       for (e = 0; e < coneSize; e++) {
1131ac17e89SToby Isaac         PetscInt point = cone[e];
1141ac17e89SToby Isaac         PetscInt refpoint = refCone[e];
1151ac17e89SToby Isaac         PetscDualSpace espace;
1161ac17e89SToby Isaac 
1179566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetPointSubspace(cspace,refpoint,&espace));
1181ac17e89SToby Isaac         if (sp->pointSpaces[point-pStart] == NULL) {
1199566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)espace));
1201ac17e89SToby Isaac           sp->pointSpaces[point-pStart] = espace;
1211ac17e89SToby Isaac         }
1221ac17e89SToby Isaac       }
1231ac17e89SToby Isaac     }
1241ac17e89SToby Isaac   }
1259566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetSection(sp, &section));
1269566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(sp, &spdim));
1279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(spdim, &(sp->functional)));
1289566063dSJacob Faibussowitsch   PetscCall(PetscDualSpacePushForwardSubspaces_Internal(sp, pStart, pEnd));
1291ac17e89SToby Isaac   PetscFunctionReturn(0);
1301ac17e89SToby Isaac }
1311ac17e89SToby Isaac 
1321ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceRefinedView_Ascii(PetscDualSpace sp, PetscViewer viewer)
1331ac17e89SToby Isaac {
1341ac17e89SToby Isaac   PetscFunctionBegin;
1351ac17e89SToby Isaac   if (sp->dm && sp->pointSpaces) {
1361ac17e89SToby Isaac     PetscInt pStart, pEnd;
1371ac17e89SToby Isaac     PetscInt cStart, cEnd, c;
1381ac17e89SToby Isaac 
1399566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(sp->dm, &pStart, &pEnd));
1409566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(sp->dm, 0, &cStart, &cEnd));
1419566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Refined dual space:\n"));
1429566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1431ac17e89SToby Isaac     for (c = cStart; c < cEnd; c++) {
1441ac17e89SToby Isaac       if (!sp->pointSpaces[c-pStart]) {
1459566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Cell space %D not set yet\n", c));
1461ac17e89SToby Isaac       } else {
1479566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Cell space %D:ot set yet\n", c));
1489566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
1499566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceView(sp->pointSpaces[c-pStart],viewer));
1509566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
1511ac17e89SToby Isaac       }
1521ac17e89SToby Isaac     }
1539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1541ac17e89SToby Isaac   } else {
1559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Refined dual space: (cell spaces not set yet)\n"));
1561ac17e89SToby Isaac   }
1571ac17e89SToby Isaac   PetscFunctionReturn(0);
1581ac17e89SToby Isaac }
1591ac17e89SToby Isaac 
1601ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceView_Refined(PetscDualSpace sp, PetscViewer viewer)
1611ac17e89SToby Isaac {
1621ac17e89SToby Isaac   PetscBool      iascii;
1631ac17e89SToby Isaac 
1641ac17e89SToby Isaac   PetscFunctionBegin;
1651ac17e89SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1661ac17e89SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1679566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
1689566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscDualSpaceRefinedView_Ascii(sp, viewer));
1691ac17e89SToby Isaac   PetscFunctionReturn(0);
1701ac17e89SToby Isaac }
1711ac17e89SToby Isaac 
1721ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceInitialize_Refined(PetscDualSpace sp)
1731ac17e89SToby Isaac {
1741ac17e89SToby Isaac   PetscFunctionBegin;
1751ac17e89SToby Isaac   sp->ops->destroy              = PetscDualSpaceDestroy_Refined;
1761ac17e89SToby Isaac   sp->ops->view                 = PetscDualSpaceView_Refined;
1771ac17e89SToby Isaac   sp->ops->setfromoptions       = NULL;
1781ac17e89SToby Isaac   sp->ops->duplicate            = NULL;
1791ac17e89SToby Isaac   sp->ops->setup                = PetscDualSpaceSetUp_Refined;
1801ac17e89SToby Isaac   sp->ops->createheightsubspace = NULL;
1811ac17e89SToby Isaac   sp->ops->createpointsubspace  = NULL;
1821ac17e89SToby Isaac   sp->ops->getsymmetries        = NULL;
1831ac17e89SToby Isaac   sp->ops->apply                = PetscDualSpaceApplyDefault;
1841ac17e89SToby Isaac   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
1851ac17e89SToby Isaac   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
1861ac17e89SToby Isaac   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
1871ac17e89SToby Isaac   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
1881ac17e89SToby Isaac   PetscFunctionReturn(0);
1891ac17e89SToby Isaac }
1901ac17e89SToby Isaac 
1911ac17e89SToby Isaac /*MC
1921ac17e89SToby Isaac   PETSCDUALSPACEREFINED = "refined" - A PetscDualSpace object that defines the joint dual space of a group of cells, usually refined from one larger cell
1931ac17e89SToby Isaac 
1941ac17e89SToby Isaac   Level: intermediate
1951ac17e89SToby Isaac 
1961ac17e89SToby Isaac .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType()
1971ac17e89SToby Isaac M*/
1981ac17e89SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Refined(PetscDualSpace sp)
1991ac17e89SToby Isaac {
2001ac17e89SToby Isaac   PetscDualSpace_Refined *ref;
2011ac17e89SToby Isaac 
2021ac17e89SToby Isaac   PetscFunctionBegin;
2031ac17e89SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2049566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(sp,&ref));
2051ac17e89SToby Isaac   sp->data = ref;
2061ac17e89SToby Isaac 
2079566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceInitialize_Refined(sp));
2089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceRefinedSetCellSpaces_C", PetscDualSpaceRefinedSetCellSpaces_Refined));
2091ac17e89SToby Isaac   PetscFunctionReturn(0);
2101ac17e89SToby Isaac }
211