xref: /petsc/src/dm/dt/dualspace/impls/refined/dualspacerefined.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode ierr;
261ac17e89SToby Isaac 
271ac17e89SToby Isaac   PetscFunctionBegin;
281ac17e89SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
291ac17e89SToby Isaac   PetscValidPointer(cellSpaces,2);
30*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change cell spaces after setup is called");
311ac17e89SToby Isaac   ierr = PetscTryMethod(sp, "PetscDualSpaceRefinedSetCellSpaces_C", (PetscDualSpace,const PetscDualSpace []),(sp,cellSpaces));CHKERRQ(ierr);
321ac17e89SToby Isaac   PetscFunctionReturn(0);
331ac17e89SToby Isaac }
341ac17e89SToby Isaac 
351ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceRefinedSetCellSpaces_Refined(PetscDualSpace sp, const PetscDualSpace cellSpaces[])
361ac17e89SToby Isaac {
371ac17e89SToby Isaac   DM dm;
381ac17e89SToby Isaac   PetscInt pStart, pEnd;
391ac17e89SToby Isaac   PetscInt cStart, cEnd, c;
401ac17e89SToby Isaac   PetscErrorCode ierr;
411ac17e89SToby Isaac 
421ac17e89SToby Isaac   PetscFunctionBegin;
431ac17e89SToby Isaac   dm = sp->dm;
44*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!dm,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_WRONGSTATE, "PetscDualSpace must have a DM (PetscDualSpaceSetDM()) before calling PetscDualSpaceRefinedSetCellSpaces");
451ac17e89SToby Isaac   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
461ac17e89SToby Isaac   if (!sp->pointSpaces) {
471ac17e89SToby Isaac 
481ac17e89SToby Isaac     ierr = PetscCalloc1(pEnd-pStart,&(sp->pointSpaces));CHKERRQ(ierr);
491ac17e89SToby Isaac   }
501ac17e89SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
511ac17e89SToby Isaac   for (c = 0; c < cEnd - cStart; c++) {
521ac17e89SToby Isaac     ierr = PetscObjectReference((PetscObject)cellSpaces[c]);CHKERRQ(ierr);
531ac17e89SToby Isaac     ierr = PetscDualSpaceDestroy(&(sp->pointSpaces[c + cStart - pStart]));CHKERRQ(ierr);
541ac17e89SToby Isaac     sp->pointSpaces[c+cStart-pStart] = cellSpaces[c];
551ac17e89SToby Isaac   }
561ac17e89SToby Isaac   PetscFunctionReturn(0);
571ac17e89SToby Isaac }
581ac17e89SToby Isaac 
591ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceDestroy_Refined(PetscDualSpace sp)
601ac17e89SToby Isaac {
611ac17e89SToby Isaac   PetscDualSpace_Refined *ref = (PetscDualSpace_Refined *) sp->data;
621ac17e89SToby Isaac   PetscErrorCode          ierr;
631ac17e89SToby Isaac 
641ac17e89SToby Isaac   PetscFunctionBegin;
651ac17e89SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceRefinedSetCellSpaces_C", NULL);CHKERRQ(ierr);
661ac17e89SToby Isaac   ierr = PetscFree(ref);CHKERRQ(ierr);
671ac17e89SToby Isaac   PetscFunctionReturn(0);
681ac17e89SToby Isaac }
691ac17e89SToby Isaac 
701ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceSetUp_Refined(PetscDualSpace sp)
711ac17e89SToby Isaac {
721ac17e89SToby Isaac   PetscInt pStart, pEnd, depth;
731ac17e89SToby Isaac   PetscInt cStart, cEnd, c, spdim;
741ac17e89SToby Isaac   PetscInt h;
751ac17e89SToby Isaac   DM dm;
761ac17e89SToby Isaac   PetscSection   section;
771ac17e89SToby Isaac   PetscErrorCode ierr;
781ac17e89SToby Isaac 
791ac17e89SToby Isaac   PetscFunctionBegin;
801ac17e89SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
811ac17e89SToby Isaac   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
821ac17e89SToby Isaac   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
831ac17e89SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
841ac17e89SToby Isaac   for (c = cStart; c < cEnd; c++) {
851ac17e89SToby Isaac     if (sp->pointSpaces[c-pStart]) {
861ac17e89SToby Isaac       PetscInt ccStart, ccEnd;
87*2c71b3e2SJacob 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");
88*2c71b3e2SJacob 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");
891ac17e89SToby Isaac       ierr = DMPlexGetHeightStratum(sp->pointSpaces[c-pStart]->dm, 0, &ccStart, &ccEnd);CHKERRQ(ierr);
90*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(ccEnd - ccStart != 1,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "All cell spaces must have a single cell themselves");
911ac17e89SToby Isaac     }
921ac17e89SToby Isaac   }
931ac17e89SToby Isaac   for (c = cStart; c < cEnd; c++) {
941ac17e89SToby Isaac     if (sp->pointSpaces[c-pStart]) {
951ac17e89SToby Isaac       PetscBool cUniform;
961ac17e89SToby Isaac 
971ac17e89SToby Isaac       ierr = PetscDualSpaceGetUniform(sp->pointSpaces[c-pStart],&cUniform);CHKERRQ(ierr);
981ac17e89SToby Isaac       if (!cUniform) break;
991ac17e89SToby Isaac     }
1001ac17e89SToby Isaac     if ((c > cStart) && sp->pointSpaces[c-pStart] != sp->pointSpaces[c-1-pStart]) break;
1011ac17e89SToby Isaac   }
1021ac17e89SToby Isaac   if (c < cEnd) sp->uniform = PETSC_FALSE;
1031ac17e89SToby Isaac   for (h = 0; h < depth; h++) {
1041ac17e89SToby Isaac     PetscInt hStart, hEnd;
1051ac17e89SToby Isaac 
1061ac17e89SToby Isaac     ierr = DMPlexGetHeightStratum(dm, h, &hStart, &hEnd);CHKERRQ(ierr);
1071ac17e89SToby Isaac     for (c = hStart; c < hEnd; c++) {
1081ac17e89SToby Isaac       PetscInt coneSize, e;
1091ac17e89SToby Isaac       PetscDualSpace cspace = sp->pointSpaces[c-pStart];
1101ac17e89SToby Isaac       const PetscInt *cone;
1111ac17e89SToby Isaac       const PetscInt *refCone;
1121ac17e89SToby Isaac 
1131ac17e89SToby Isaac       if (!cspace) continue;
1141ac17e89SToby Isaac       ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
1151ac17e89SToby Isaac       ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
1161ac17e89SToby Isaac       ierr = DMPlexGetCone(cspace->dm, 0, &refCone);CHKERRQ(ierr);
1171ac17e89SToby Isaac       for (e = 0; e < coneSize; e++) {
1181ac17e89SToby Isaac         PetscInt point = cone[e];
1191ac17e89SToby Isaac         PetscInt refpoint = refCone[e];
1201ac17e89SToby Isaac         PetscDualSpace espace;
1211ac17e89SToby Isaac 
1221ac17e89SToby Isaac         ierr = PetscDualSpaceGetPointSubspace(cspace,refpoint,&espace);CHKERRQ(ierr);
1231ac17e89SToby Isaac         if (sp->pointSpaces[point-pStart] == NULL) {
1241ac17e89SToby Isaac           ierr = PetscObjectReference((PetscObject)espace);CHKERRQ(ierr);
1251ac17e89SToby Isaac           sp->pointSpaces[point-pStart] = espace;
1261ac17e89SToby Isaac         }
1271ac17e89SToby Isaac       }
1281ac17e89SToby Isaac     }
1291ac17e89SToby Isaac   }
1301ac17e89SToby Isaac   ierr = PetscDualSpaceGetSection(sp, &section);CHKERRQ(ierr);
1311ac17e89SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr);
1321ac17e89SToby Isaac   ierr = PetscMalloc1(spdim, &(sp->functional));CHKERRQ(ierr);
1331ac17e89SToby Isaac   ierr = PetscDualSpacePushForwardSubspaces_Internal(sp, pStart, pEnd);CHKERRQ(ierr);
1341ac17e89SToby Isaac   PetscFunctionReturn(0);
1351ac17e89SToby Isaac }
1361ac17e89SToby Isaac 
1371ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceRefinedView_Ascii(PetscDualSpace sp, PetscViewer viewer)
1381ac17e89SToby Isaac {
1391ac17e89SToby Isaac   PetscErrorCode      ierr;
1401ac17e89SToby Isaac 
1411ac17e89SToby Isaac   PetscFunctionBegin;
1421ac17e89SToby Isaac   if (sp->dm && sp->pointSpaces) {
1431ac17e89SToby Isaac     PetscInt pStart, pEnd;
1441ac17e89SToby Isaac     PetscInt cStart, cEnd, c;
1451ac17e89SToby Isaac 
1461ac17e89SToby Isaac     ierr = DMPlexGetChart(sp->dm, &pStart, &pEnd);CHKERRQ(ierr);
1471ac17e89SToby Isaac     ierr = DMPlexGetHeightStratum(sp->dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1481ac17e89SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "Refined dual space:\n");CHKERRQ(ierr);
1491ac17e89SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1501ac17e89SToby Isaac     for (c = cStart; c < cEnd; c++) {
1511ac17e89SToby Isaac       if (!sp->pointSpaces[c-pStart]) {
1521ac17e89SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Cell space %D not set yet\n", c);CHKERRQ(ierr);
1531ac17e89SToby Isaac       } else {
1541ac17e89SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Cell space %D:ot set yet\n", c);CHKERRQ(ierr);
1551ac17e89SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1561ac17e89SToby Isaac         ierr = PetscDualSpaceView(sp->pointSpaces[c-pStart],viewer);CHKERRQ(ierr);
1571ac17e89SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1581ac17e89SToby Isaac       }
1591ac17e89SToby Isaac     }
1601ac17e89SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1611ac17e89SToby Isaac   } else {
1621ac17e89SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "Refined dual space: (cell spaces not set yet)\n");CHKERRQ(ierr);
1631ac17e89SToby Isaac   }
1641ac17e89SToby Isaac   PetscFunctionReturn(0);
1651ac17e89SToby Isaac }
1661ac17e89SToby Isaac 
1671ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceView_Refined(PetscDualSpace sp, PetscViewer viewer)
1681ac17e89SToby Isaac {
1691ac17e89SToby Isaac   PetscBool      iascii;
1701ac17e89SToby Isaac   PetscErrorCode ierr;
1711ac17e89SToby Isaac 
1721ac17e89SToby Isaac   PetscFunctionBegin;
1731ac17e89SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1741ac17e89SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1751ac17e89SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
1761ac17e89SToby Isaac   if (iascii) {ierr = PetscDualSpaceRefinedView_Ascii(sp, viewer);CHKERRQ(ierr);}
1771ac17e89SToby Isaac   PetscFunctionReturn(0);
1781ac17e89SToby Isaac }
1791ac17e89SToby Isaac 
1801ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceInitialize_Refined(PetscDualSpace sp)
1811ac17e89SToby Isaac {
1821ac17e89SToby Isaac   PetscFunctionBegin;
1831ac17e89SToby Isaac   sp->ops->destroy              = PetscDualSpaceDestroy_Refined;
1841ac17e89SToby Isaac   sp->ops->view                 = PetscDualSpaceView_Refined;
1851ac17e89SToby Isaac   sp->ops->setfromoptions       = NULL;
1861ac17e89SToby Isaac   sp->ops->duplicate            = NULL;
1871ac17e89SToby Isaac   sp->ops->setup                = PetscDualSpaceSetUp_Refined;
1881ac17e89SToby Isaac   sp->ops->createheightsubspace = NULL;
1891ac17e89SToby Isaac   sp->ops->createpointsubspace  = NULL;
1901ac17e89SToby Isaac   sp->ops->getsymmetries        = NULL;
1911ac17e89SToby Isaac   sp->ops->apply                = PetscDualSpaceApplyDefault;
1921ac17e89SToby Isaac   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
1931ac17e89SToby Isaac   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
1941ac17e89SToby Isaac   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
1951ac17e89SToby Isaac   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
1961ac17e89SToby Isaac   PetscFunctionReturn(0);
1971ac17e89SToby Isaac }
1981ac17e89SToby Isaac 
1991ac17e89SToby Isaac /*MC
2001ac17e89SToby Isaac   PETSCDUALSPACEREFINED = "refined" - A PetscDualSpace object that defines the joint dual space of a group of cells, usually refined from one larger cell
2011ac17e89SToby Isaac 
2021ac17e89SToby Isaac   Level: intermediate
2031ac17e89SToby Isaac 
2041ac17e89SToby Isaac .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType()
2051ac17e89SToby Isaac M*/
2061ac17e89SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Refined(PetscDualSpace sp)
2071ac17e89SToby Isaac {
2081ac17e89SToby Isaac   PetscDualSpace_Refined *ref;
2091ac17e89SToby Isaac   PetscErrorCode      ierr;
2101ac17e89SToby Isaac 
2111ac17e89SToby Isaac   PetscFunctionBegin;
2121ac17e89SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2131ac17e89SToby Isaac   ierr     = PetscNewLog(sp,&ref);CHKERRQ(ierr);
2141ac17e89SToby Isaac   sp->data = ref;
2151ac17e89SToby Isaac 
2161ac17e89SToby Isaac   ierr = PetscDualSpaceInitialize_Refined(sp);CHKERRQ(ierr);
2171ac17e89SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceRefinedSetCellSpaces_C", PetscDualSpaceRefinedSetCellSpaces_Refined);CHKERRQ(ierr);
2181ac17e89SToby Isaac   PetscFunctionReturn(0);
2191ac17e89SToby Isaac }
220