1*1ac17e89SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2*1ac17e89SToby Isaac #include <petscdmplex.h> 3*1ac17e89SToby Isaac 4*1ac17e89SToby Isaac typedef struct { 5*1ac17e89SToby Isaac PetscInt dummy; 6*1ac17e89SToby Isaac } PetscDualSpace_Refined; 7*1ac17e89SToby Isaac 8*1ac17e89SToby Isaac /*@ 9*1ac17e89SToby Isaac PetscDualSpaceRefinedSetCellSpaces - Set the dual spaces for the closures of each of the cells 10*1ac17e89SToby Isaac in the multicell DM of a PetscDualSpace 11*1ac17e89SToby Isaac 12*1ac17e89SToby Isaac Collective on PetscDualSpace 13*1ac17e89SToby Isaac 14*1ac17e89SToby Isaac Input Arguments: 15*1ac17e89SToby Isaac + sp - a PetscDualSpace 16*1ac17e89SToby Isaac - cellSpaces - one PetscDualSpace for each of the cells. The reference count of each cell space will be incremented, 17*1ac17e89SToby Isaac so the user is still responsible for these spaces afterwards 18*1ac17e89SToby Isaac 19*1ac17e89SToby Isaac Level: intermediate 20*1ac17e89SToby Isaac 21*1ac17e89SToby Isaac .seealso: PetscFERefine() 22*1ac17e89SToby Isaac @*/ 23*1ac17e89SToby Isaac PetscErrorCode PetscDualSpaceRefinedSetCellSpaces(PetscDualSpace sp, const PetscDualSpace cellSpaces[]) 24*1ac17e89SToby Isaac { 25*1ac17e89SToby Isaac PetscErrorCode ierr; 26*1ac17e89SToby Isaac 27*1ac17e89SToby Isaac PetscFunctionBegin; 28*1ac17e89SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 29*1ac17e89SToby Isaac PetscValidPointer(cellSpaces,2); 30*1ac17e89SToby Isaac if (sp->setupcalled) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change cell spaces after setup is called"); 31*1ac17e89SToby Isaac ierr = PetscTryMethod(sp, "PetscDualSpaceRefinedSetCellSpaces_C", (PetscDualSpace,const PetscDualSpace []),(sp,cellSpaces));CHKERRQ(ierr); 32*1ac17e89SToby Isaac PetscFunctionReturn(0); 33*1ac17e89SToby Isaac } 34*1ac17e89SToby Isaac 35*1ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceRefinedSetCellSpaces_Refined(PetscDualSpace sp, const PetscDualSpace cellSpaces[]) 36*1ac17e89SToby Isaac { 37*1ac17e89SToby Isaac DM dm; 38*1ac17e89SToby Isaac PetscInt pStart, pEnd; 39*1ac17e89SToby Isaac PetscInt cStart, cEnd, c; 40*1ac17e89SToby Isaac PetscErrorCode ierr; 41*1ac17e89SToby Isaac 42*1ac17e89SToby Isaac PetscFunctionBegin; 43*1ac17e89SToby Isaac dm = sp->dm; 44*1ac17e89SToby Isaac if (!dm) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_WRONGSTATE, "PetscDualSpace must have a DM (PetscDualSpaceSetDM()) before calling PetscDualSpaceRefinedSetCellSpaces"); 45*1ac17e89SToby Isaac ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 46*1ac17e89SToby Isaac if (!sp->pointSpaces) { 47*1ac17e89SToby Isaac 48*1ac17e89SToby Isaac ierr = PetscCalloc1(pEnd-pStart,&(sp->pointSpaces));CHKERRQ(ierr); 49*1ac17e89SToby Isaac } 50*1ac17e89SToby Isaac ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 51*1ac17e89SToby Isaac for (c = 0; c < cEnd - cStart; c++) { 52*1ac17e89SToby Isaac ierr = PetscObjectReference((PetscObject)cellSpaces[c]);CHKERRQ(ierr); 53*1ac17e89SToby Isaac ierr = PetscDualSpaceDestroy(&(sp->pointSpaces[c + cStart - pStart]));CHKERRQ(ierr); 54*1ac17e89SToby Isaac sp->pointSpaces[c+cStart-pStart] = cellSpaces[c]; 55*1ac17e89SToby Isaac } 56*1ac17e89SToby Isaac PetscFunctionReturn(0); 57*1ac17e89SToby Isaac } 58*1ac17e89SToby Isaac 59*1ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceDestroy_Refined(PetscDualSpace sp) 60*1ac17e89SToby Isaac { 61*1ac17e89SToby Isaac PetscDualSpace_Refined *ref = (PetscDualSpace_Refined *) sp->data; 62*1ac17e89SToby Isaac PetscErrorCode ierr; 63*1ac17e89SToby Isaac 64*1ac17e89SToby Isaac PetscFunctionBegin; 65*1ac17e89SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceRefinedSetCellSpaces_C", NULL);CHKERRQ(ierr); 66*1ac17e89SToby Isaac ierr = PetscFree(ref);CHKERRQ(ierr); 67*1ac17e89SToby Isaac PetscFunctionReturn(0); 68*1ac17e89SToby Isaac } 69*1ac17e89SToby Isaac 70*1ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceSetUp_Refined(PetscDualSpace sp) 71*1ac17e89SToby Isaac { 72*1ac17e89SToby Isaac PetscInt pStart, pEnd, depth; 73*1ac17e89SToby Isaac PetscInt cStart, cEnd, c, spdim; 74*1ac17e89SToby Isaac PetscInt h; 75*1ac17e89SToby Isaac DM dm; 76*1ac17e89SToby Isaac PetscSection section; 77*1ac17e89SToby Isaac PetscErrorCode ierr; 78*1ac17e89SToby Isaac 79*1ac17e89SToby Isaac PetscFunctionBegin; 80*1ac17e89SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 81*1ac17e89SToby Isaac ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 82*1ac17e89SToby Isaac ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 83*1ac17e89SToby Isaac ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 84*1ac17e89SToby Isaac for (c = cStart; c < cEnd; c++) { 85*1ac17e89SToby Isaac if (sp->pointSpaces[c-pStart]) { 86*1ac17e89SToby Isaac PetscInt ccStart, ccEnd; 87*1ac17e89SToby Isaac if (sp->pointSpaces[c-pStart]->k != sp->k) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "All cell spaces must have the same form degree as the refined dual space"); 88*1ac17e89SToby Isaac if (sp->pointSpaces[c-pStart]->Nc != sp->Nc) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "All cell spaces must have the same number of components as the refined dual space"); 89*1ac17e89SToby Isaac ierr = DMPlexGetHeightStratum(sp->pointSpaces[c-pStart]->dm, 0, &ccStart, &ccEnd);CHKERRQ(ierr); 90*1ac17e89SToby Isaac if (ccEnd - ccStart != 1) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "All cell spaces must have a single cell themselves"); 91*1ac17e89SToby Isaac } 92*1ac17e89SToby Isaac } 93*1ac17e89SToby Isaac for (c = cStart; c < cEnd; c++) { 94*1ac17e89SToby Isaac if (sp->pointSpaces[c-pStart]) { 95*1ac17e89SToby Isaac PetscBool cUniform; 96*1ac17e89SToby Isaac 97*1ac17e89SToby Isaac ierr = PetscDualSpaceGetUniform(sp->pointSpaces[c-pStart],&cUniform);CHKERRQ(ierr); 98*1ac17e89SToby Isaac if (!cUniform) break; 99*1ac17e89SToby Isaac } 100*1ac17e89SToby Isaac if ((c > cStart) && sp->pointSpaces[c-pStart] != sp->pointSpaces[c-1-pStart]) break; 101*1ac17e89SToby Isaac } 102*1ac17e89SToby Isaac if (c < cEnd) sp->uniform = PETSC_FALSE; 103*1ac17e89SToby Isaac for (h = 0; h < depth; h++) { 104*1ac17e89SToby Isaac PetscInt hStart, hEnd; 105*1ac17e89SToby Isaac 106*1ac17e89SToby Isaac ierr = DMPlexGetHeightStratum(dm, h, &hStart, &hEnd);CHKERRQ(ierr); 107*1ac17e89SToby Isaac for (c = hStart; c < hEnd; c++) { 108*1ac17e89SToby Isaac PetscInt coneSize, e; 109*1ac17e89SToby Isaac PetscDualSpace cspace = sp->pointSpaces[c-pStart]; 110*1ac17e89SToby Isaac const PetscInt *cone; 111*1ac17e89SToby Isaac const PetscInt *refCone; 112*1ac17e89SToby Isaac 113*1ac17e89SToby Isaac if (!cspace) continue; 114*1ac17e89SToby Isaac ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 115*1ac17e89SToby Isaac ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 116*1ac17e89SToby Isaac ierr = DMPlexGetCone(cspace->dm, 0, &refCone);CHKERRQ(ierr); 117*1ac17e89SToby Isaac for (e = 0; e < coneSize; e++) { 118*1ac17e89SToby Isaac PetscInt point = cone[e]; 119*1ac17e89SToby Isaac PetscInt refpoint = refCone[e]; 120*1ac17e89SToby Isaac PetscDualSpace espace; 121*1ac17e89SToby Isaac 122*1ac17e89SToby Isaac ierr = PetscDualSpaceGetPointSubspace(cspace,refpoint,&espace);CHKERRQ(ierr); 123*1ac17e89SToby Isaac if (sp->pointSpaces[point-pStart] == NULL) { 124*1ac17e89SToby Isaac ierr = PetscObjectReference((PetscObject)espace);CHKERRQ(ierr); 125*1ac17e89SToby Isaac sp->pointSpaces[point-pStart] = espace; 126*1ac17e89SToby Isaac } 127*1ac17e89SToby Isaac } 128*1ac17e89SToby Isaac } 129*1ac17e89SToby Isaac } 130*1ac17e89SToby Isaac ierr = PetscDualSpaceGetSection(sp, §ion);CHKERRQ(ierr); 131*1ac17e89SToby Isaac ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr); 132*1ac17e89SToby Isaac ierr = PetscMalloc1(spdim, &(sp->functional));CHKERRQ(ierr); 133*1ac17e89SToby Isaac ierr = PetscDualSpacePushForwardSubspaces_Internal(sp, pStart, pEnd);CHKERRQ(ierr); 134*1ac17e89SToby Isaac PetscFunctionReturn(0); 135*1ac17e89SToby Isaac } 136*1ac17e89SToby Isaac 137*1ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceRefinedView_Ascii(PetscDualSpace sp, PetscViewer viewer) 138*1ac17e89SToby Isaac { 139*1ac17e89SToby Isaac PetscErrorCode ierr; 140*1ac17e89SToby Isaac 141*1ac17e89SToby Isaac PetscFunctionBegin; 142*1ac17e89SToby Isaac if (sp->dm && sp->pointSpaces) { 143*1ac17e89SToby Isaac PetscInt pStart, pEnd; 144*1ac17e89SToby Isaac PetscInt cStart, cEnd, c; 145*1ac17e89SToby Isaac 146*1ac17e89SToby Isaac ierr = DMPlexGetChart(sp->dm, &pStart, &pEnd);CHKERRQ(ierr); 147*1ac17e89SToby Isaac ierr = DMPlexGetHeightStratum(sp->dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 148*1ac17e89SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Refined dual space:\n");CHKERRQ(ierr); 149*1ac17e89SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 150*1ac17e89SToby Isaac for (c = cStart; c < cEnd; c++) { 151*1ac17e89SToby Isaac if (!sp->pointSpaces[c-pStart]) { 152*1ac17e89SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Cell space %D not set yet\n", c);CHKERRQ(ierr); 153*1ac17e89SToby Isaac } else { 154*1ac17e89SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Cell space %D:ot set yet\n", c);CHKERRQ(ierr); 155*1ac17e89SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 156*1ac17e89SToby Isaac ierr = PetscDualSpaceView(sp->pointSpaces[c-pStart],viewer);CHKERRQ(ierr); 157*1ac17e89SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 158*1ac17e89SToby Isaac } 159*1ac17e89SToby Isaac } 160*1ac17e89SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 161*1ac17e89SToby Isaac } else { 162*1ac17e89SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Refined dual space: (cell spaces not set yet)\n");CHKERRQ(ierr); 163*1ac17e89SToby Isaac } 164*1ac17e89SToby Isaac PetscFunctionReturn(0); 165*1ac17e89SToby Isaac } 166*1ac17e89SToby Isaac 167*1ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceView_Refined(PetscDualSpace sp, PetscViewer viewer) 168*1ac17e89SToby Isaac { 169*1ac17e89SToby Isaac PetscBool iascii; 170*1ac17e89SToby Isaac PetscErrorCode ierr; 171*1ac17e89SToby Isaac 172*1ac17e89SToby Isaac PetscFunctionBegin; 173*1ac17e89SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 174*1ac17e89SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 175*1ac17e89SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 176*1ac17e89SToby Isaac if (iascii) {ierr = PetscDualSpaceRefinedView_Ascii(sp, viewer);CHKERRQ(ierr);} 177*1ac17e89SToby Isaac PetscFunctionReturn(0); 178*1ac17e89SToby Isaac } 179*1ac17e89SToby Isaac 180*1ac17e89SToby Isaac static PetscErrorCode PetscDualSpaceInitialize_Refined(PetscDualSpace sp) 181*1ac17e89SToby Isaac { 182*1ac17e89SToby Isaac PetscFunctionBegin; 183*1ac17e89SToby Isaac sp->ops->destroy = PetscDualSpaceDestroy_Refined; 184*1ac17e89SToby Isaac sp->ops->view = PetscDualSpaceView_Refined; 185*1ac17e89SToby Isaac sp->ops->setfromoptions = NULL; 186*1ac17e89SToby Isaac sp->ops->duplicate = NULL; 187*1ac17e89SToby Isaac sp->ops->setup = PetscDualSpaceSetUp_Refined; 188*1ac17e89SToby Isaac sp->ops->createheightsubspace = NULL; 189*1ac17e89SToby Isaac sp->ops->createpointsubspace = NULL; 190*1ac17e89SToby Isaac sp->ops->getsymmetries = NULL; 191*1ac17e89SToby Isaac sp->ops->apply = PetscDualSpaceApplyDefault; 192*1ac17e89SToby Isaac sp->ops->applyall = PetscDualSpaceApplyAllDefault; 193*1ac17e89SToby Isaac sp->ops->applyint = PetscDualSpaceApplyInteriorDefault; 194*1ac17e89SToby Isaac sp->ops->createalldata = PetscDualSpaceCreateAllDataDefault; 195*1ac17e89SToby Isaac sp->ops->createintdata = PetscDualSpaceCreateInteriorDataDefault; 196*1ac17e89SToby Isaac PetscFunctionReturn(0); 197*1ac17e89SToby Isaac } 198*1ac17e89SToby Isaac 199*1ac17e89SToby Isaac /*MC 200*1ac17e89SToby Isaac PETSCDUALSPACEREFINED = "refined" - A PetscDualSpace object that defines the joint dual space of a group of cells, usually refined from one larger cell 201*1ac17e89SToby Isaac 202*1ac17e89SToby Isaac Level: intermediate 203*1ac17e89SToby Isaac 204*1ac17e89SToby Isaac .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 205*1ac17e89SToby Isaac M*/ 206*1ac17e89SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Refined(PetscDualSpace sp) 207*1ac17e89SToby Isaac { 208*1ac17e89SToby Isaac PetscDualSpace_Refined *ref; 209*1ac17e89SToby Isaac PetscErrorCode ierr; 210*1ac17e89SToby Isaac 211*1ac17e89SToby Isaac PetscFunctionBegin; 212*1ac17e89SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 213*1ac17e89SToby Isaac ierr = PetscNewLog(sp,&ref);CHKERRQ(ierr); 214*1ac17e89SToby Isaac sp->data = ref; 215*1ac17e89SToby Isaac 216*1ac17e89SToby Isaac ierr = PetscDualSpaceInitialize_Refined(sp);CHKERRQ(ierr); 217*1ac17e89SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceRefinedSetCellSpaces_C", PetscDualSpaceRefinedSetCellSpaces_Refined);CHKERRQ(ierr); 218*1ac17e89SToby Isaac PetscFunctionReturn(0); 219*1ac17e89SToby Isaac } 220