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 21db781477SPatrick Sanan .seealso: `PetscFERefine()` 221ac17e89SToby Isaac @*/ 239371c9d4SSatish Balay PetscErrorCode PetscDualSpaceRefinedSetCellSpaces(PetscDualSpace sp, const PetscDualSpace cellSpaces[]) { 241ac17e89SToby Isaac PetscFunctionBegin; 251ac17e89SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 261ac17e89SToby Isaac PetscValidPointer(cellSpaces, 2); 2728b400f6SJacob Faibussowitsch PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change cell spaces after setup is called"); 28cac4c232SBarry Smith PetscTryMethod(sp, "PetscDualSpaceRefinedSetCellSpaces_C", (PetscDualSpace, const PetscDualSpace[]), (sp, cellSpaces)); 291ac17e89SToby Isaac PetscFunctionReturn(0); 301ac17e89SToby Isaac } 311ac17e89SToby Isaac 329371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceRefinedSetCellSpaces_Refined(PetscDualSpace sp, const PetscDualSpace cellSpaces[]) { 331ac17e89SToby Isaac DM dm; 341ac17e89SToby Isaac PetscInt pStart, pEnd; 351ac17e89SToby Isaac PetscInt cStart, cEnd, c; 361ac17e89SToby Isaac 371ac17e89SToby Isaac PetscFunctionBegin; 381ac17e89SToby Isaac dm = sp->dm; 3928b400f6SJacob Faibussowitsch PetscCheck(dm, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "PetscDualSpace must have a DM (PetscDualSpaceSetDM()) before calling PetscDualSpaceRefinedSetCellSpaces"); 409566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 4148a46eb9SPierre Jolivet if (!sp->pointSpaces) PetscCall(PetscCalloc1(pEnd - pStart, &(sp->pointSpaces))); 429566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 431ac17e89SToby Isaac for (c = 0; c < cEnd - cStart; c++) { 449566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)cellSpaces[c])); 459566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&(sp->pointSpaces[c + cStart - pStart]))); 461ac17e89SToby Isaac sp->pointSpaces[c + cStart - pStart] = cellSpaces[c]; 471ac17e89SToby Isaac } 481ac17e89SToby Isaac PetscFunctionReturn(0); 491ac17e89SToby Isaac } 501ac17e89SToby Isaac 519371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceDestroy_Refined(PetscDualSpace sp) { 521ac17e89SToby Isaac PetscDualSpace_Refined *ref = (PetscDualSpace_Refined *)sp->data; 531ac17e89SToby Isaac 541ac17e89SToby Isaac PetscFunctionBegin; 559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceRefinedSetCellSpaces_C", NULL)); 569566063dSJacob Faibussowitsch PetscCall(PetscFree(ref)); 571ac17e89SToby Isaac PetscFunctionReturn(0); 581ac17e89SToby Isaac } 591ac17e89SToby Isaac 609371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceSetUp_Refined(PetscDualSpace sp) { 611ac17e89SToby Isaac PetscInt pStart, pEnd, depth; 621ac17e89SToby Isaac PetscInt cStart, cEnd, c, spdim; 631ac17e89SToby Isaac PetscInt h; 641ac17e89SToby Isaac DM dm; 651ac17e89SToby Isaac PetscSection section; 661ac17e89SToby Isaac 671ac17e89SToby Isaac PetscFunctionBegin; 689566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp, &dm)); 699566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 709566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 719566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 721ac17e89SToby Isaac for (c = cStart; c < cEnd; c++) { 731ac17e89SToby Isaac if (sp->pointSpaces[c - pStart]) { 741ac17e89SToby Isaac PetscInt ccStart, ccEnd; 7508401ef6SPierre 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"); 7608401ef6SPierre 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"); 779566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(sp->pointSpaces[c - pStart]->dm, 0, &ccStart, &ccEnd)); 781dca8a05SBarry Smith PetscCheck(ccEnd - ccStart == 1, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "All cell spaces must have a single cell themselves"); 791ac17e89SToby Isaac } 801ac17e89SToby Isaac } 811ac17e89SToby Isaac for (c = cStart; c < cEnd; c++) { 821ac17e89SToby Isaac if (sp->pointSpaces[c - pStart]) { 831ac17e89SToby Isaac PetscBool cUniform; 841ac17e89SToby Isaac 859566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetUniform(sp->pointSpaces[c - pStart], &cUniform)); 861ac17e89SToby Isaac if (!cUniform) break; 871ac17e89SToby Isaac } 881ac17e89SToby Isaac if ((c > cStart) && sp->pointSpaces[c - pStart] != sp->pointSpaces[c - 1 - pStart]) break; 891ac17e89SToby Isaac } 901ac17e89SToby Isaac if (c < cEnd) sp->uniform = PETSC_FALSE; 911ac17e89SToby Isaac for (h = 0; h < depth; h++) { 921ac17e89SToby Isaac PetscInt hStart, hEnd; 931ac17e89SToby Isaac 949566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd)); 951ac17e89SToby Isaac for (c = hStart; c < hEnd; c++) { 961ac17e89SToby Isaac PetscInt coneSize, e; 971ac17e89SToby Isaac PetscDualSpace cspace = sp->pointSpaces[c - pStart]; 981ac17e89SToby Isaac const PetscInt *cone; 991ac17e89SToby Isaac const PetscInt *refCone; 1001ac17e89SToby Isaac 1011ac17e89SToby Isaac if (!cspace) continue; 1029566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, c, &coneSize)); 1039566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, c, &cone)); 1049566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(cspace->dm, 0, &refCone)); 1051ac17e89SToby Isaac for (e = 0; e < coneSize; e++) { 1061ac17e89SToby Isaac PetscInt point = cone[e]; 1071ac17e89SToby Isaac PetscInt refpoint = refCone[e]; 1081ac17e89SToby Isaac PetscDualSpace espace; 1091ac17e89SToby Isaac 1109566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetPointSubspace(cspace, refpoint, &espace)); 1111ac17e89SToby Isaac if (sp->pointSpaces[point - pStart] == NULL) { 1129566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)espace)); 1131ac17e89SToby Isaac sp->pointSpaces[point - pStart] = espace; 1141ac17e89SToby Isaac } 1151ac17e89SToby Isaac } 1161ac17e89SToby Isaac } 1171ac17e89SToby Isaac } 1189566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetSection(sp, §ion)); 1199566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp, &spdim)); 1209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(spdim, &(sp->functional))); 1219566063dSJacob Faibussowitsch PetscCall(PetscDualSpacePushForwardSubspaces_Internal(sp, pStart, pEnd)); 1221ac17e89SToby Isaac PetscFunctionReturn(0); 1231ac17e89SToby Isaac } 1241ac17e89SToby Isaac 1259371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceRefinedView_Ascii(PetscDualSpace sp, PetscViewer viewer) { 1261ac17e89SToby Isaac PetscFunctionBegin; 1271ac17e89SToby Isaac if (sp->dm && sp->pointSpaces) { 1281ac17e89SToby Isaac PetscInt pStart, pEnd; 1291ac17e89SToby Isaac PetscInt cStart, cEnd, c; 1301ac17e89SToby Isaac 1319566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(sp->dm, &pStart, &pEnd)); 1329566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(sp->dm, 0, &cStart, &cEnd)); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Refined dual space:\n")); 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1351ac17e89SToby Isaac for (c = cStart; c < cEnd; c++) { 1361ac17e89SToby Isaac if (!sp->pointSpaces[c - pStart]) { 13763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Cell space %" PetscInt_FMT " not set yet\n", c)); 1381ac17e89SToby Isaac } else { 13963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Cell space %" PetscInt_FMT ":ot set yet\n", c)); 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1419566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceView(sp->pointSpaces[c - pStart], viewer)); 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1431ac17e89SToby Isaac } 1441ac17e89SToby Isaac } 1459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1461baa6e33SBarry Smith } else PetscCall(PetscViewerASCIIPrintf(viewer, "Refined dual space: (cell spaces not set yet)\n")); 1471ac17e89SToby Isaac PetscFunctionReturn(0); 1481ac17e89SToby Isaac } 1491ac17e89SToby Isaac 1509371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceView_Refined(PetscDualSpace sp, PetscViewer viewer) { 1511ac17e89SToby Isaac PetscBool iascii; 1521ac17e89SToby Isaac 1531ac17e89SToby Isaac PetscFunctionBegin; 1541ac17e89SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1551ac17e89SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1569566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1579566063dSJacob Faibussowitsch if (iascii) PetscCall(PetscDualSpaceRefinedView_Ascii(sp, viewer)); 1581ac17e89SToby Isaac PetscFunctionReturn(0); 1591ac17e89SToby Isaac } 1601ac17e89SToby Isaac 1619371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceInitialize_Refined(PetscDualSpace sp) { 1621ac17e89SToby Isaac PetscFunctionBegin; 1631ac17e89SToby Isaac sp->ops->destroy = PetscDualSpaceDestroy_Refined; 1641ac17e89SToby Isaac sp->ops->view = PetscDualSpaceView_Refined; 1651ac17e89SToby Isaac sp->ops->setfromoptions = NULL; 1661ac17e89SToby Isaac sp->ops->duplicate = NULL; 1671ac17e89SToby Isaac sp->ops->setup = PetscDualSpaceSetUp_Refined; 1681ac17e89SToby Isaac sp->ops->createheightsubspace = NULL; 1691ac17e89SToby Isaac sp->ops->createpointsubspace = NULL; 1701ac17e89SToby Isaac sp->ops->getsymmetries = NULL; 1711ac17e89SToby Isaac sp->ops->apply = PetscDualSpaceApplyDefault; 1721ac17e89SToby Isaac sp->ops->applyall = PetscDualSpaceApplyAllDefault; 1731ac17e89SToby Isaac sp->ops->applyint = PetscDualSpaceApplyInteriorDefault; 1741ac17e89SToby Isaac sp->ops->createalldata = PetscDualSpaceCreateAllDataDefault; 1751ac17e89SToby Isaac sp->ops->createintdata = PetscDualSpaceCreateInteriorDataDefault; 1761ac17e89SToby Isaac PetscFunctionReturn(0); 1771ac17e89SToby Isaac } 1781ac17e89SToby Isaac 1791ac17e89SToby Isaac /*MC 1801ac17e89SToby Isaac PETSCDUALSPACEREFINED = "refined" - A PetscDualSpace object that defines the joint dual space of a group of cells, usually refined from one larger cell 1811ac17e89SToby Isaac 1821ac17e89SToby Isaac Level: intermediate 1831ac17e89SToby Isaac 184db781477SPatrick Sanan .seealso: `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()` 1851ac17e89SToby Isaac M*/ 1869371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Refined(PetscDualSpace sp) { 1871ac17e89SToby Isaac PetscDualSpace_Refined *ref; 1881ac17e89SToby Isaac 1891ac17e89SToby Isaac PetscFunctionBegin; 1901ac17e89SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 191*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ref)); 1921ac17e89SToby Isaac sp->data = ref; 1931ac17e89SToby Isaac 1949566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceInitialize_Refined(sp)); 1959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceRefinedSetCellSpaces_C", PetscDualSpaceRefinedSetCellSpaces_Refined)); 1961ac17e89SToby Isaac PetscFunctionReturn(0); 1971ac17e89SToby Isaac } 198