120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac #include <petscdmplex.h> 320cf1dd8SToby Isaac 420cf1dd8SToby Isaac PetscClassId PETSCDUALSPACE_CLASSID = 0; 520cf1dd8SToby Isaac 6ead873ccSMatthew G. Knepley PetscLogEvent PETSCDUALSPACE_SetUp; 7ead873ccSMatthew G. Knepley 820cf1dd8SToby Isaac PetscFunctionList PetscDualSpaceList = NULL; 920cf1dd8SToby Isaac PetscBool PetscDualSpaceRegisterAllCalled = PETSC_FALSE; 1020cf1dd8SToby Isaac 116f905325SMatthew G. Knepley /* 126f905325SMatthew G. Knepley PetscDualSpaceLatticePointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to at most 'max'. 136f905325SMatthew G. Knepley Ordering is lexicographic with lowest index as least significant in ordering. 14b4457527SToby Isaac e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {0,2}. 156f905325SMatthew G. Knepley 166f905325SMatthew G. Knepley Input Parameters: 176f905325SMatthew G. Knepley + len - The length of the tuple 186f905325SMatthew G. Knepley . max - The maximum sum 196f905325SMatthew G. Knepley - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition 206f905325SMatthew G. Knepley 216f905325SMatthew G. Knepley Output Parameter: 226f905325SMatthew G. Knepley . tup - A tuple of len integers whos sum is at most 'max' 236f905325SMatthew G. Knepley 246f905325SMatthew G. Knepley Level: developer 256f905325SMatthew G. Knepley 266f905325SMatthew G. Knepley .seealso: PetscDualSpaceTensorPointLexicographic_Internal() 276f905325SMatthew G. Knepley */ 286f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceLatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[]) 296f905325SMatthew G. Knepley { 306f905325SMatthew G. Knepley PetscFunctionBegin; 316f905325SMatthew G. Knepley while (len--) { 326f905325SMatthew G. Knepley max -= tup[len]; 336f905325SMatthew G. Knepley if (!max) { 346f905325SMatthew G. Knepley tup[len] = 0; 356f905325SMatthew G. Knepley break; 366f905325SMatthew G. Knepley } 376f905325SMatthew G. Knepley } 386f905325SMatthew G. Knepley tup[++len]++; 396f905325SMatthew G. Knepley PetscFunctionReturn(0); 406f905325SMatthew G. Knepley } 416f905325SMatthew G. Knepley 426f905325SMatthew G. Knepley /* 436f905325SMatthew G. Knepley PetscDualSpaceTensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'. 446f905325SMatthew G. Knepley Ordering is lexicographic with lowest index as least significant in ordering. 456f905325SMatthew G. Knepley e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}. 466f905325SMatthew G. Knepley 476f905325SMatthew G. Knepley Input Parameters: 486f905325SMatthew G. Knepley + len - The length of the tuple 496f905325SMatthew G. Knepley . max - The maximum value 506f905325SMatthew G. Knepley - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition 516f905325SMatthew G. Knepley 526f905325SMatthew G. Knepley Output Parameter: 536f905325SMatthew G. Knepley . tup - A tuple of len integers whos sum is at most 'max' 546f905325SMatthew G. Knepley 556f905325SMatthew G. Knepley Level: developer 566f905325SMatthew G. Knepley 576f905325SMatthew G. Knepley .seealso: PetscDualSpaceLatticePointLexicographic_Internal() 586f905325SMatthew G. Knepley */ 596f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceTensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[]) 606f905325SMatthew G. Knepley { 616f905325SMatthew G. Knepley PetscInt i; 626f905325SMatthew G. Knepley 636f905325SMatthew G. Knepley PetscFunctionBegin; 646f905325SMatthew G. Knepley for (i = 0; i < len; i++) { 656f905325SMatthew G. Knepley if (tup[i] < max) { 666f905325SMatthew G. Knepley break; 676f905325SMatthew G. Knepley } else { 686f905325SMatthew G. Knepley tup[i] = 0; 696f905325SMatthew G. Knepley } 706f905325SMatthew G. Knepley } 716f905325SMatthew G. Knepley tup[i]++; 726f905325SMatthew G. Knepley PetscFunctionReturn(0); 736f905325SMatthew G. Knepley } 746f905325SMatthew G. Knepley 7520cf1dd8SToby Isaac /*@C 7620cf1dd8SToby Isaac PetscDualSpaceRegister - Adds a new PetscDualSpace implementation 7720cf1dd8SToby Isaac 7820cf1dd8SToby Isaac Not Collective 7920cf1dd8SToby Isaac 8020cf1dd8SToby Isaac Input Parameters: 8120cf1dd8SToby Isaac + name - The name of a new user-defined creation routine 8220cf1dd8SToby Isaac - create_func - The creation routine itself 8320cf1dd8SToby Isaac 8420cf1dd8SToby Isaac Notes: 8520cf1dd8SToby Isaac PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces 8620cf1dd8SToby Isaac 8720cf1dd8SToby Isaac Sample usage: 8820cf1dd8SToby Isaac .vb 8920cf1dd8SToby Isaac PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate); 9020cf1dd8SToby Isaac .ve 9120cf1dd8SToby Isaac 9220cf1dd8SToby Isaac Then, your PetscDualSpace type can be chosen with the procedural interface via 9320cf1dd8SToby Isaac .vb 9420cf1dd8SToby Isaac PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *); 9520cf1dd8SToby Isaac PetscDualSpaceSetType(PetscDualSpace, "my_dual_space"); 9620cf1dd8SToby Isaac .ve 9720cf1dd8SToby Isaac or at runtime via the option 9820cf1dd8SToby Isaac .vb 9920cf1dd8SToby Isaac -petscdualspace_type my_dual_space 10020cf1dd8SToby Isaac .ve 10120cf1dd8SToby Isaac 10220cf1dd8SToby Isaac Level: advanced 10320cf1dd8SToby Isaac 10420cf1dd8SToby Isaac .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy() 10520cf1dd8SToby Isaac 10620cf1dd8SToby Isaac @*/ 10720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace)) 10820cf1dd8SToby Isaac { 10920cf1dd8SToby Isaac PetscErrorCode ierr; 11020cf1dd8SToby Isaac 11120cf1dd8SToby Isaac PetscFunctionBegin; 11220cf1dd8SToby Isaac ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr); 11320cf1dd8SToby Isaac PetscFunctionReturn(0); 11420cf1dd8SToby Isaac } 11520cf1dd8SToby Isaac 11620cf1dd8SToby Isaac /*@C 11720cf1dd8SToby Isaac PetscDualSpaceSetType - Builds a particular PetscDualSpace 11820cf1dd8SToby Isaac 119d083f849SBarry Smith Collective on sp 12020cf1dd8SToby Isaac 12120cf1dd8SToby Isaac Input Parameters: 12220cf1dd8SToby Isaac + sp - The PetscDualSpace object 12320cf1dd8SToby Isaac - name - The kind of space 12420cf1dd8SToby Isaac 12520cf1dd8SToby Isaac Options Database Key: 12620cf1dd8SToby Isaac . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types 12720cf1dd8SToby Isaac 12820cf1dd8SToby Isaac Level: intermediate 12920cf1dd8SToby Isaac 13020cf1dd8SToby Isaac .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate() 13120cf1dd8SToby Isaac @*/ 13220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name) 13320cf1dd8SToby Isaac { 13420cf1dd8SToby Isaac PetscErrorCode (*r)(PetscDualSpace); 13520cf1dd8SToby Isaac PetscBool match; 13620cf1dd8SToby Isaac PetscErrorCode ierr; 13720cf1dd8SToby Isaac 13820cf1dd8SToby Isaac PetscFunctionBegin; 13920cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 14020cf1dd8SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr); 14120cf1dd8SToby Isaac if (match) PetscFunctionReturn(0); 14220cf1dd8SToby Isaac 14320cf1dd8SToby Isaac if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);} 14420cf1dd8SToby Isaac ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr); 1452c71b3e2SJacob Faibussowitsch PetscCheckFalse(!r,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name); 14620cf1dd8SToby Isaac 14720cf1dd8SToby Isaac if (sp->ops->destroy) { 14820cf1dd8SToby Isaac ierr = (*sp->ops->destroy)(sp);CHKERRQ(ierr); 14920cf1dd8SToby Isaac sp->ops->destroy = NULL; 15020cf1dd8SToby Isaac } 15120cf1dd8SToby Isaac ierr = (*r)(sp);CHKERRQ(ierr); 15220cf1dd8SToby Isaac ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr); 15320cf1dd8SToby Isaac PetscFunctionReturn(0); 15420cf1dd8SToby Isaac } 15520cf1dd8SToby Isaac 15620cf1dd8SToby Isaac /*@C 15720cf1dd8SToby Isaac PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object. 15820cf1dd8SToby Isaac 15920cf1dd8SToby Isaac Not Collective 16020cf1dd8SToby Isaac 16120cf1dd8SToby Isaac Input Parameter: 16220cf1dd8SToby Isaac . sp - The PetscDualSpace 16320cf1dd8SToby Isaac 16420cf1dd8SToby Isaac Output Parameter: 16520cf1dd8SToby Isaac . name - The PetscDualSpace type name 16620cf1dd8SToby Isaac 16720cf1dd8SToby Isaac Level: intermediate 16820cf1dd8SToby Isaac 16920cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate() 17020cf1dd8SToby Isaac @*/ 17120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name) 17220cf1dd8SToby Isaac { 17320cf1dd8SToby Isaac PetscErrorCode ierr; 17420cf1dd8SToby Isaac 17520cf1dd8SToby Isaac PetscFunctionBegin; 17620cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 17720cf1dd8SToby Isaac PetscValidPointer(name, 2); 17820cf1dd8SToby Isaac if (!PetscDualSpaceRegisterAllCalled) { 17920cf1dd8SToby Isaac ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr); 18020cf1dd8SToby Isaac } 18120cf1dd8SToby Isaac *name = ((PetscObject) sp)->type_name; 18220cf1dd8SToby Isaac PetscFunctionReturn(0); 18320cf1dd8SToby Isaac } 18420cf1dd8SToby Isaac 185221d6281SMatthew G. Knepley static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v) 186221d6281SMatthew G. Knepley { 187221d6281SMatthew G. Knepley PetscViewerFormat format; 188221d6281SMatthew G. Knepley PetscInt pdim, f; 189221d6281SMatthew G. Knepley PetscErrorCode ierr; 190221d6281SMatthew G. Knepley 191221d6281SMatthew G. Knepley PetscFunctionBegin; 192221d6281SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 193221d6281SMatthew G. Knepley ierr = PetscObjectPrintClassNamePrefixType((PetscObject) sp, v);CHKERRQ(ierr); 194221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 195b4457527SToby Isaac if (sp->k) { 196b4457527SToby Isaac ierr = PetscViewerASCIIPrintf(v, "Dual space for %D-forms %swith %D components, size %D\n", PetscAbsInt(sp->k), sp->k < 0 ? "(stored in dual form) ": "", sp->Nc, pdim);CHKERRQ(ierr); 197b4457527SToby Isaac } else { 198221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim);CHKERRQ(ierr); 199b4457527SToby Isaac } 200221d6281SMatthew G. Knepley if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);} 201221d6281SMatthew G. Knepley ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr); 202221d6281SMatthew G. Knepley if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 203221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 204221d6281SMatthew G. Knepley for (f = 0; f < pdim; ++f) { 205221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f);CHKERRQ(ierr); 206221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 207221d6281SMatthew G. Knepley ierr = PetscQuadratureView(sp->functional[f], v);CHKERRQ(ierr); 208221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 209221d6281SMatthew G. Knepley } 210221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 211221d6281SMatthew G. Knepley } 212221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 213221d6281SMatthew G. Knepley PetscFunctionReturn(0); 214221d6281SMatthew G. Knepley } 215221d6281SMatthew G. Knepley 216fe2efc57SMark /*@C 217fe2efc57SMark PetscDualSpaceViewFromOptions - View from Options 218fe2efc57SMark 219fe2efc57SMark Collective on PetscDualSpace 220fe2efc57SMark 221fe2efc57SMark Input Parameters: 222fe2efc57SMark + A - the PetscDualSpace object 223736c3998SJose E. Roman . obj - Optional object, proivides prefix 224736c3998SJose E. Roman - name - command line option 225fe2efc57SMark 226fe2efc57SMark Level: intermediate 227fe2efc57SMark .seealso: PetscDualSpace, PetscDualSpaceView(), PetscObjectViewFromOptions(), PetscDualSpaceCreate() 228fe2efc57SMark @*/ 229fe2efc57SMark PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace A,PetscObject obj,const char name[]) 230fe2efc57SMark { 231fe2efc57SMark PetscErrorCode ierr; 232fe2efc57SMark 233fe2efc57SMark PetscFunctionBegin; 234fe2efc57SMark PetscValidHeaderSpecific(A,PETSCDUALSPACE_CLASSID,1); 235fe2efc57SMark ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr); 236fe2efc57SMark PetscFunctionReturn(0); 237fe2efc57SMark } 238fe2efc57SMark 23920cf1dd8SToby Isaac /*@ 24020cf1dd8SToby Isaac PetscDualSpaceView - Views a PetscDualSpace 24120cf1dd8SToby Isaac 242d083f849SBarry Smith Collective on sp 24320cf1dd8SToby Isaac 244d8d19677SJose E. Roman Input Parameters: 24520cf1dd8SToby Isaac + sp - the PetscDualSpace object to view 24620cf1dd8SToby Isaac - v - the viewer 24720cf1dd8SToby Isaac 248a4ce7ad1SMatthew G. Knepley Level: beginner 24920cf1dd8SToby Isaac 250fe2efc57SMark .seealso PetscDualSpaceDestroy(), PetscDualSpace 25120cf1dd8SToby Isaac @*/ 25220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v) 25320cf1dd8SToby Isaac { 254d9bac1caSLisandro Dalcin PetscBool iascii; 25520cf1dd8SToby Isaac PetscErrorCode ierr; 25620cf1dd8SToby Isaac 25720cf1dd8SToby Isaac PetscFunctionBegin; 25820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 259d9bac1caSLisandro Dalcin if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2); 26020cf1dd8SToby Isaac if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);} 261d9bac1caSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 262221d6281SMatthew G. Knepley if (iascii) {ierr = PetscDualSpaceView_ASCII(sp, v);CHKERRQ(ierr);} 26320cf1dd8SToby Isaac PetscFunctionReturn(0); 26420cf1dd8SToby Isaac } 26520cf1dd8SToby Isaac 26620cf1dd8SToby Isaac /*@ 26720cf1dd8SToby Isaac PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database 26820cf1dd8SToby Isaac 269d083f849SBarry Smith Collective on sp 27020cf1dd8SToby Isaac 27120cf1dd8SToby Isaac Input Parameter: 27220cf1dd8SToby Isaac . sp - the PetscDualSpace object to set options for 27320cf1dd8SToby Isaac 27420cf1dd8SToby Isaac Options Database: 2758f2aacc6SMatthew G. Knepley + -petscdualspace_order <order> - the approximation order of the space 276fe36a153SMatthew G. Knepley . -petscdualspace_form_degree <deg> - the form degree, say 0 for point evaluations, or 2 for area integrals 2778f2aacc6SMatthew G. Knepley . -petscdualspace_components <c> - the number of components, say d for a vector field 2788f2aacc6SMatthew G. Knepley - -petscdualspace_refcell <celltype> - Reference cell type name 27920cf1dd8SToby Isaac 280a4ce7ad1SMatthew G. Knepley Level: intermediate 28120cf1dd8SToby Isaac 282fe2efc57SMark .seealso PetscDualSpaceView(), PetscDualSpace, PetscObjectSetFromOptions() 28320cf1dd8SToby Isaac @*/ 28420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp) 28520cf1dd8SToby Isaac { 2862df84da0SMatthew G. Knepley DMPolytopeType refCell = DM_POLYTOPE_TRIANGLE; 28720cf1dd8SToby Isaac const char *defaultType; 28820cf1dd8SToby Isaac char name[256]; 289f783ec47SMatthew G. Knepley PetscBool flg; 29020cf1dd8SToby Isaac PetscErrorCode ierr; 29120cf1dd8SToby Isaac 29220cf1dd8SToby Isaac PetscFunctionBegin; 29320cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 29420cf1dd8SToby Isaac if (!((PetscObject) sp)->type_name) { 29520cf1dd8SToby Isaac defaultType = PETSCDUALSPACELAGRANGE; 29620cf1dd8SToby Isaac } else { 29720cf1dd8SToby Isaac defaultType = ((PetscObject) sp)->type_name; 29820cf1dd8SToby Isaac } 29920cf1dd8SToby Isaac if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 30020cf1dd8SToby Isaac 30120cf1dd8SToby Isaac ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 30220cf1dd8SToby Isaac ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr); 30320cf1dd8SToby Isaac if (flg) { 30420cf1dd8SToby Isaac ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr); 30520cf1dd8SToby Isaac } else if (!((PetscObject) sp)->type_name) { 30620cf1dd8SToby Isaac ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr); 30720cf1dd8SToby Isaac } 308b4457527SToby Isaac ierr = PetscOptionsBoundedInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL,0);CHKERRQ(ierr); 309b4457527SToby Isaac ierr = PetscOptionsInt("-petscdualspace_form_degree", "The form degree of the dofs", "PetscDualSpaceSetFormDegree", sp->k, &sp->k, NULL);CHKERRQ(ierr); 3105a856986SBarry Smith ierr = PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL,1);CHKERRQ(ierr); 31120cf1dd8SToby Isaac if (sp->ops->setfromoptions) { 31220cf1dd8SToby Isaac ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr); 31320cf1dd8SToby Isaac } 314f783ec47SMatthew G. Knepley ierr = PetscOptionsEnum("-petscdualspace_refcell", "Reference cell shape", "PetscDualSpaceSetReferenceCell", DMPolytopeTypes, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);CHKERRQ(ierr); 315063ee4adSMatthew G. Knepley if (flg) { 316063ee4adSMatthew G. Knepley DM K; 317063ee4adSMatthew G. Knepley 318f783ec47SMatthew G. Knepley ierr = DMPlexCreateReferenceCell(PETSC_COMM_SELF, refCell, &K);CHKERRQ(ierr); 319063ee4adSMatthew G. Knepley ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr); 320063ee4adSMatthew G. Knepley ierr = DMDestroy(&K);CHKERRQ(ierr); 321063ee4adSMatthew G. Knepley } 322063ee4adSMatthew G. Knepley 32320cf1dd8SToby Isaac /* process any options handlers added with PetscObjectAddOptionsHandler() */ 32420cf1dd8SToby Isaac ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr); 32520cf1dd8SToby Isaac ierr = PetscOptionsEnd();CHKERRQ(ierr); 326063ee4adSMatthew G. Knepley sp->setfromoptionscalled = PETSC_TRUE; 32720cf1dd8SToby Isaac PetscFunctionReturn(0); 32820cf1dd8SToby Isaac } 32920cf1dd8SToby Isaac 33020cf1dd8SToby Isaac /*@ 33120cf1dd8SToby Isaac PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace 33220cf1dd8SToby Isaac 333d083f849SBarry Smith Collective on sp 33420cf1dd8SToby Isaac 33520cf1dd8SToby Isaac Input Parameter: 33620cf1dd8SToby Isaac . sp - the PetscDualSpace object to setup 33720cf1dd8SToby Isaac 338a4ce7ad1SMatthew G. Knepley Level: intermediate 33920cf1dd8SToby Isaac 340fe2efc57SMark .seealso PetscDualSpaceView(), PetscDualSpaceDestroy(), PetscDualSpace 34120cf1dd8SToby Isaac @*/ 34220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp) 34320cf1dd8SToby Isaac { 34420cf1dd8SToby Isaac PetscErrorCode ierr; 34520cf1dd8SToby Isaac 34620cf1dd8SToby Isaac PetscFunctionBegin; 34720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 34820cf1dd8SToby Isaac if (sp->setupcalled) PetscFunctionReturn(0); 349ead873ccSMatthew G. Knepley ierr = PetscLogEventBegin(PETSCDUALSPACE_SetUp, sp, 0, 0, 0);CHKERRQ(ierr); 35020cf1dd8SToby Isaac sp->setupcalled = PETSC_TRUE; 35120cf1dd8SToby Isaac if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);} 352ead873ccSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCDUALSPACE_SetUp, sp, 0, 0, 0);CHKERRQ(ierr); 353063ee4adSMatthew G. Knepley if (sp->setfromoptionscalled) {ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);} 35420cf1dd8SToby Isaac PetscFunctionReturn(0); 35520cf1dd8SToby Isaac } 35620cf1dd8SToby Isaac 357b4457527SToby Isaac static PetscErrorCode PetscDualSpaceClearDMData_Internal(PetscDualSpace sp, DM dm) 358b4457527SToby Isaac { 359b4457527SToby Isaac PetscInt pStart = -1, pEnd = -1, depth = -1; 360b4457527SToby Isaac PetscErrorCode ierr; 361b4457527SToby Isaac 362b4457527SToby Isaac PetscFunctionBegin; 363b4457527SToby Isaac if (!dm) PetscFunctionReturn(0); 364b4457527SToby Isaac ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 365b4457527SToby Isaac ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 366b4457527SToby Isaac 367b4457527SToby Isaac if (sp->pointSpaces) { 368b4457527SToby Isaac PetscInt i; 369b4457527SToby Isaac 370b4457527SToby Isaac for (i = 0; i < pEnd - pStart; i++) { 371b4457527SToby Isaac ierr = PetscDualSpaceDestroy(&(sp->pointSpaces[i]));CHKERRQ(ierr); 372b4457527SToby Isaac } 373b4457527SToby Isaac } 374b4457527SToby Isaac ierr = PetscFree(sp->pointSpaces);CHKERRQ(ierr); 375b4457527SToby Isaac 376b4457527SToby Isaac if (sp->heightSpaces) { 377b4457527SToby Isaac PetscInt i; 378b4457527SToby Isaac 379b4457527SToby Isaac for (i = 0; i <= depth; i++) { 380b4457527SToby Isaac ierr = PetscDualSpaceDestroy(&(sp->heightSpaces[i]));CHKERRQ(ierr); 381b4457527SToby Isaac } 382b4457527SToby Isaac } 383b4457527SToby Isaac ierr = PetscFree(sp->heightSpaces);CHKERRQ(ierr); 384b4457527SToby Isaac 385b4457527SToby Isaac ierr = PetscSectionDestroy(&(sp->pointSection));CHKERRQ(ierr); 386b4457527SToby Isaac ierr = PetscQuadratureDestroy(&(sp->intNodes));CHKERRQ(ierr); 387b4457527SToby Isaac ierr = VecDestroy(&(sp->intDofValues));CHKERRQ(ierr); 388b4457527SToby Isaac ierr = VecDestroy(&(sp->intNodeValues));CHKERRQ(ierr); 389b4457527SToby Isaac ierr = MatDestroy(&(sp->intMat));CHKERRQ(ierr); 390b4457527SToby Isaac ierr = PetscQuadratureDestroy(&(sp->allNodes));CHKERRQ(ierr); 391b4457527SToby Isaac ierr = VecDestroy(&(sp->allDofValues));CHKERRQ(ierr); 392b4457527SToby Isaac ierr = VecDestroy(&(sp->allNodeValues));CHKERRQ(ierr); 393b4457527SToby Isaac ierr = MatDestroy(&(sp->allMat));CHKERRQ(ierr); 394b4457527SToby Isaac ierr = PetscFree(sp->numDof);CHKERRQ(ierr); 395b4457527SToby Isaac PetscFunctionReturn(0); 396b4457527SToby Isaac } 397b4457527SToby Isaac 39820cf1dd8SToby Isaac /*@ 39920cf1dd8SToby Isaac PetscDualSpaceDestroy - Destroys a PetscDualSpace object 40020cf1dd8SToby Isaac 401d083f849SBarry Smith Collective on sp 40220cf1dd8SToby Isaac 40320cf1dd8SToby Isaac Input Parameter: 40420cf1dd8SToby Isaac . sp - the PetscDualSpace object to destroy 40520cf1dd8SToby Isaac 406a4ce7ad1SMatthew G. Knepley Level: beginner 40720cf1dd8SToby Isaac 408fe2efc57SMark .seealso PetscDualSpaceView(), PetscDualSpace(), PetscDualSpaceCreate() 40920cf1dd8SToby Isaac @*/ 41020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp) 41120cf1dd8SToby Isaac { 41220cf1dd8SToby Isaac PetscInt dim, f; 413b4457527SToby Isaac DM dm; 41420cf1dd8SToby Isaac PetscErrorCode ierr; 41520cf1dd8SToby Isaac 41620cf1dd8SToby Isaac PetscFunctionBegin; 41720cf1dd8SToby Isaac if (!*sp) PetscFunctionReturn(0); 41820cf1dd8SToby Isaac PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1); 41920cf1dd8SToby Isaac 420ea78f98cSLisandro Dalcin if (--((PetscObject)(*sp))->refct > 0) {*sp = NULL; PetscFunctionReturn(0);} 42120cf1dd8SToby Isaac ((PetscObject) (*sp))->refct = 0; 42220cf1dd8SToby Isaac 42320cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr); 424b4457527SToby Isaac dm = (*sp)->dm; 425b4457527SToby Isaac 426b4457527SToby Isaac if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);} 427b4457527SToby Isaac ierr = PetscDualSpaceClearDMData_Internal(*sp, dm);CHKERRQ(ierr); 428b4457527SToby Isaac 42920cf1dd8SToby Isaac for (f = 0; f < dim; ++f) { 43020cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr); 43120cf1dd8SToby Isaac } 43220cf1dd8SToby Isaac ierr = PetscFree((*sp)->functional);CHKERRQ(ierr); 43320cf1dd8SToby Isaac ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr); 43420cf1dd8SToby Isaac ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr); 43520cf1dd8SToby Isaac PetscFunctionReturn(0); 43620cf1dd8SToby Isaac } 43720cf1dd8SToby Isaac 43820cf1dd8SToby Isaac /*@ 43920cf1dd8SToby Isaac PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType(). 44020cf1dd8SToby Isaac 441d083f849SBarry Smith Collective 44220cf1dd8SToby Isaac 44320cf1dd8SToby Isaac Input Parameter: 44420cf1dd8SToby Isaac . comm - The communicator for the PetscDualSpace object 44520cf1dd8SToby Isaac 44620cf1dd8SToby Isaac Output Parameter: 44720cf1dd8SToby Isaac . sp - The PetscDualSpace object 44820cf1dd8SToby Isaac 44920cf1dd8SToby Isaac Level: beginner 45020cf1dd8SToby Isaac 45120cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE 45220cf1dd8SToby Isaac @*/ 45320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp) 45420cf1dd8SToby Isaac { 45520cf1dd8SToby Isaac PetscDualSpace s; 45620cf1dd8SToby Isaac PetscErrorCode ierr; 45720cf1dd8SToby Isaac 45820cf1dd8SToby Isaac PetscFunctionBegin; 45920cf1dd8SToby Isaac PetscValidPointer(sp, 2); 46020cf1dd8SToby Isaac ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr); 46120cf1dd8SToby Isaac *sp = NULL; 46220cf1dd8SToby Isaac ierr = PetscFEInitializePackage();CHKERRQ(ierr); 46320cf1dd8SToby Isaac 46420cf1dd8SToby Isaac ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr); 46520cf1dd8SToby Isaac 46620cf1dd8SToby Isaac s->order = 0; 46720cf1dd8SToby Isaac s->Nc = 1; 4684bee2e38SMatthew G. Knepley s->k = 0; 469b4457527SToby Isaac s->spdim = -1; 470b4457527SToby Isaac s->spintdim = -1; 471b4457527SToby Isaac s->uniform = PETSC_TRUE; 47220cf1dd8SToby Isaac s->setupcalled = PETSC_FALSE; 47320cf1dd8SToby Isaac 47420cf1dd8SToby Isaac *sp = s; 47520cf1dd8SToby Isaac PetscFunctionReturn(0); 47620cf1dd8SToby Isaac } 47720cf1dd8SToby Isaac 47820cf1dd8SToby Isaac /*@ 47920cf1dd8SToby Isaac PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup. 48020cf1dd8SToby Isaac 481d083f849SBarry Smith Collective on sp 48220cf1dd8SToby Isaac 48320cf1dd8SToby Isaac Input Parameter: 48420cf1dd8SToby Isaac . sp - The original PetscDualSpace 48520cf1dd8SToby Isaac 48620cf1dd8SToby Isaac Output Parameter: 48720cf1dd8SToby Isaac . spNew - The duplicate PetscDualSpace 48820cf1dd8SToby Isaac 48920cf1dd8SToby Isaac Level: beginner 49020cf1dd8SToby Isaac 49120cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType() 49220cf1dd8SToby Isaac @*/ 49320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew) 49420cf1dd8SToby Isaac { 495b4457527SToby Isaac DM dm; 496b4457527SToby Isaac PetscDualSpaceType type; 497b4457527SToby Isaac const char *name; 49820cf1dd8SToby Isaac PetscErrorCode ierr; 49920cf1dd8SToby Isaac 50020cf1dd8SToby Isaac PetscFunctionBegin; 50120cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 50220cf1dd8SToby Isaac PetscValidPointer(spNew, 2); 503b4457527SToby Isaac ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject)sp), spNew);CHKERRQ(ierr); 504b4457527SToby Isaac ierr = PetscObjectGetName((PetscObject) sp, &name);CHKERRQ(ierr); 505b4457527SToby Isaac ierr = PetscObjectSetName((PetscObject) *spNew, name);CHKERRQ(ierr); 506b4457527SToby Isaac ierr = PetscDualSpaceGetType(sp, &type);CHKERRQ(ierr); 507b4457527SToby Isaac ierr = PetscDualSpaceSetType(*spNew, type);CHKERRQ(ierr); 508b4457527SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 509b4457527SToby Isaac ierr = PetscDualSpaceSetDM(*spNew, dm);CHKERRQ(ierr); 510b4457527SToby Isaac 511b4457527SToby Isaac (*spNew)->order = sp->order; 512b4457527SToby Isaac (*spNew)->k = sp->k; 513b4457527SToby Isaac (*spNew)->Nc = sp->Nc; 514b4457527SToby Isaac (*spNew)->uniform = sp->uniform; 515b4457527SToby Isaac if (sp->ops->duplicate) {ierr = (*sp->ops->duplicate)(sp, *spNew);CHKERRQ(ierr);} 51620cf1dd8SToby Isaac PetscFunctionReturn(0); 51720cf1dd8SToby Isaac } 51820cf1dd8SToby Isaac 51920cf1dd8SToby Isaac /*@ 52020cf1dd8SToby Isaac PetscDualSpaceGetDM - Get the DM representing the reference cell 52120cf1dd8SToby Isaac 52220cf1dd8SToby Isaac Not collective 52320cf1dd8SToby Isaac 52420cf1dd8SToby Isaac Input Parameter: 52520cf1dd8SToby Isaac . sp - The PetscDualSpace 52620cf1dd8SToby Isaac 52720cf1dd8SToby Isaac Output Parameter: 52820cf1dd8SToby Isaac . dm - The reference cell 52920cf1dd8SToby Isaac 53020cf1dd8SToby Isaac Level: intermediate 53120cf1dd8SToby Isaac 53220cf1dd8SToby Isaac .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate() 53320cf1dd8SToby Isaac @*/ 53420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm) 53520cf1dd8SToby Isaac { 53620cf1dd8SToby Isaac PetscFunctionBegin; 53720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 53820cf1dd8SToby Isaac PetscValidPointer(dm, 2); 53920cf1dd8SToby Isaac *dm = sp->dm; 54020cf1dd8SToby Isaac PetscFunctionReturn(0); 54120cf1dd8SToby Isaac } 54220cf1dd8SToby Isaac 54320cf1dd8SToby Isaac /*@ 54420cf1dd8SToby Isaac PetscDualSpaceSetDM - Get the DM representing the reference cell 54520cf1dd8SToby Isaac 54620cf1dd8SToby Isaac Not collective 54720cf1dd8SToby Isaac 54820cf1dd8SToby Isaac Input Parameters: 54920cf1dd8SToby Isaac + sp - The PetscDualSpace 55020cf1dd8SToby Isaac - dm - The reference cell 55120cf1dd8SToby Isaac 55220cf1dd8SToby Isaac Level: intermediate 55320cf1dd8SToby Isaac 55420cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate() 55520cf1dd8SToby Isaac @*/ 55620cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm) 55720cf1dd8SToby Isaac { 55820cf1dd8SToby Isaac PetscErrorCode ierr; 55920cf1dd8SToby Isaac 56020cf1dd8SToby Isaac PetscFunctionBegin; 56120cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 56220cf1dd8SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 5632c71b3e2SJacob Faibussowitsch PetscCheckFalse(sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change DM after dualspace is set up"); 56420cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 565b4457527SToby Isaac if (sp->dm && sp->dm != dm) { 566b4457527SToby Isaac ierr = PetscDualSpaceClearDMData_Internal(sp, sp->dm);CHKERRQ(ierr); 567b4457527SToby Isaac } 568b4457527SToby Isaac ierr = DMDestroy(&sp->dm);CHKERRQ(ierr); 56920cf1dd8SToby Isaac sp->dm = dm; 57020cf1dd8SToby Isaac PetscFunctionReturn(0); 57120cf1dd8SToby Isaac } 57220cf1dd8SToby Isaac 57320cf1dd8SToby Isaac /*@ 57420cf1dd8SToby Isaac PetscDualSpaceGetOrder - Get the order of the dual space 57520cf1dd8SToby Isaac 57620cf1dd8SToby Isaac Not collective 57720cf1dd8SToby Isaac 57820cf1dd8SToby Isaac Input Parameter: 57920cf1dd8SToby Isaac . sp - The PetscDualSpace 58020cf1dd8SToby Isaac 58120cf1dd8SToby Isaac Output Parameter: 58220cf1dd8SToby Isaac . order - The order 58320cf1dd8SToby Isaac 58420cf1dd8SToby Isaac Level: intermediate 58520cf1dd8SToby Isaac 58620cf1dd8SToby Isaac .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate() 58720cf1dd8SToby Isaac @*/ 58820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order) 58920cf1dd8SToby Isaac { 59020cf1dd8SToby Isaac PetscFunctionBegin; 59120cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 59220cf1dd8SToby Isaac PetscValidPointer(order, 2); 59320cf1dd8SToby Isaac *order = sp->order; 59420cf1dd8SToby Isaac PetscFunctionReturn(0); 59520cf1dd8SToby Isaac } 59620cf1dd8SToby Isaac 59720cf1dd8SToby Isaac /*@ 59820cf1dd8SToby Isaac PetscDualSpaceSetOrder - Set the order of the dual space 59920cf1dd8SToby Isaac 60020cf1dd8SToby Isaac Not collective 60120cf1dd8SToby Isaac 60220cf1dd8SToby Isaac Input Parameters: 60320cf1dd8SToby Isaac + sp - The PetscDualSpace 60420cf1dd8SToby Isaac - order - The order 60520cf1dd8SToby Isaac 60620cf1dd8SToby Isaac Level: intermediate 60720cf1dd8SToby Isaac 60820cf1dd8SToby Isaac .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate() 60920cf1dd8SToby Isaac @*/ 61020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order) 61120cf1dd8SToby Isaac { 61220cf1dd8SToby Isaac PetscFunctionBegin; 61320cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 6142c71b3e2SJacob Faibussowitsch PetscCheckFalse(sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change order after dualspace is set up"); 61520cf1dd8SToby Isaac sp->order = order; 61620cf1dd8SToby Isaac PetscFunctionReturn(0); 61720cf1dd8SToby Isaac } 61820cf1dd8SToby Isaac 61920cf1dd8SToby Isaac /*@ 62020cf1dd8SToby Isaac PetscDualSpaceGetNumComponents - Return the number of components for this space 62120cf1dd8SToby Isaac 62220cf1dd8SToby Isaac Input Parameter: 62320cf1dd8SToby Isaac . sp - The PetscDualSpace 62420cf1dd8SToby Isaac 62520cf1dd8SToby Isaac Output Parameter: 62620cf1dd8SToby Isaac . Nc - The number of components 62720cf1dd8SToby Isaac 62820cf1dd8SToby Isaac Note: A vector space, for example, will have d components, where d is the spatial dimension 62920cf1dd8SToby Isaac 63020cf1dd8SToby Isaac Level: intermediate 63120cf1dd8SToby Isaac 63220cf1dd8SToby Isaac .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace 63320cf1dd8SToby Isaac @*/ 63420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc) 63520cf1dd8SToby Isaac { 63620cf1dd8SToby Isaac PetscFunctionBegin; 63720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 63820cf1dd8SToby Isaac PetscValidPointer(Nc, 2); 63920cf1dd8SToby Isaac *Nc = sp->Nc; 64020cf1dd8SToby Isaac PetscFunctionReturn(0); 64120cf1dd8SToby Isaac } 64220cf1dd8SToby Isaac 64320cf1dd8SToby Isaac /*@ 64420cf1dd8SToby Isaac PetscDualSpaceSetNumComponents - Set the number of components for this space 64520cf1dd8SToby Isaac 64620cf1dd8SToby Isaac Input Parameters: 64720cf1dd8SToby Isaac + sp - The PetscDualSpace 64820cf1dd8SToby Isaac - order - The number of components 64920cf1dd8SToby Isaac 65020cf1dd8SToby Isaac Level: intermediate 65120cf1dd8SToby Isaac 65220cf1dd8SToby Isaac .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace 65320cf1dd8SToby Isaac @*/ 65420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc) 65520cf1dd8SToby Isaac { 65620cf1dd8SToby Isaac PetscFunctionBegin; 65720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 6582c71b3e2SJacob Faibussowitsch PetscCheckFalse(sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up"); 65920cf1dd8SToby Isaac sp->Nc = Nc; 66020cf1dd8SToby Isaac PetscFunctionReturn(0); 66120cf1dd8SToby Isaac } 66220cf1dd8SToby Isaac 66320cf1dd8SToby Isaac /*@ 66420cf1dd8SToby Isaac PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space 66520cf1dd8SToby Isaac 66620cf1dd8SToby Isaac Not collective 66720cf1dd8SToby Isaac 66820cf1dd8SToby Isaac Input Parameters: 66920cf1dd8SToby Isaac + sp - The PetscDualSpace 67020cf1dd8SToby Isaac - i - The basis number 67120cf1dd8SToby Isaac 67220cf1dd8SToby Isaac Output Parameter: 67320cf1dd8SToby Isaac . functional - The basis functional 67420cf1dd8SToby Isaac 67520cf1dd8SToby Isaac Level: intermediate 67620cf1dd8SToby Isaac 67720cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate() 67820cf1dd8SToby Isaac @*/ 67920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional) 68020cf1dd8SToby Isaac { 68120cf1dd8SToby Isaac PetscInt dim; 68220cf1dd8SToby Isaac PetscErrorCode ierr; 68320cf1dd8SToby Isaac 68420cf1dd8SToby Isaac PetscFunctionBegin; 68520cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 68620cf1dd8SToby Isaac PetscValidPointer(functional, 3); 68720cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr); 6882c71b3e2SJacob Faibussowitsch PetscCheckFalse((i < 0) || (i >= dim),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim); 68920cf1dd8SToby Isaac *functional = sp->functional[i]; 69020cf1dd8SToby Isaac PetscFunctionReturn(0); 69120cf1dd8SToby Isaac } 69220cf1dd8SToby Isaac 69320cf1dd8SToby Isaac /*@ 69420cf1dd8SToby Isaac PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals 69520cf1dd8SToby Isaac 69620cf1dd8SToby Isaac Not collective 69720cf1dd8SToby Isaac 69820cf1dd8SToby Isaac Input Parameter: 69920cf1dd8SToby Isaac . sp - The PetscDualSpace 70020cf1dd8SToby Isaac 70120cf1dd8SToby Isaac Output Parameter: 70220cf1dd8SToby Isaac . dim - The dimension 70320cf1dd8SToby Isaac 70420cf1dd8SToby Isaac Level: intermediate 70520cf1dd8SToby Isaac 70620cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 70720cf1dd8SToby Isaac @*/ 70820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim) 70920cf1dd8SToby Isaac { 71020cf1dd8SToby Isaac PetscErrorCode ierr; 71120cf1dd8SToby Isaac 71220cf1dd8SToby Isaac PetscFunctionBegin; 71320cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 71420cf1dd8SToby Isaac PetscValidPointer(dim, 2); 715b4457527SToby Isaac if (sp->spdim < 0) { 716b4457527SToby Isaac PetscSection section; 717b4457527SToby Isaac 718b4457527SToby Isaac ierr = PetscDualSpaceGetSection(sp, §ion);CHKERRQ(ierr); 719b4457527SToby Isaac if (section) { 720b4457527SToby Isaac ierr = PetscSectionGetStorageSize(section, &(sp->spdim));CHKERRQ(ierr); 721b4457527SToby Isaac } else sp->spdim = 0; 722b4457527SToby Isaac } 723b4457527SToby Isaac *dim = sp->spdim; 72420cf1dd8SToby Isaac PetscFunctionReturn(0); 72520cf1dd8SToby Isaac } 72620cf1dd8SToby Isaac 727b4457527SToby Isaac /*@ 728b4457527SToby Isaac PetscDualSpaceGetInteriorDimension - Get the interior dimension of the dual space, i.e. the number of basis functionals assigned to the interior of the reference domain 729b4457527SToby Isaac 730b4457527SToby Isaac Not collective 731b4457527SToby Isaac 732b4457527SToby Isaac Input Parameter: 733b4457527SToby Isaac . sp - The PetscDualSpace 734b4457527SToby Isaac 735b4457527SToby Isaac Output Parameter: 736b4457527SToby Isaac . dim - The dimension 737b4457527SToby Isaac 738b4457527SToby Isaac Level: intermediate 739b4457527SToby Isaac 740b4457527SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 741b4457527SToby Isaac @*/ 742b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace sp, PetscInt *intdim) 743b4457527SToby Isaac { 744b4457527SToby Isaac PetscErrorCode ierr; 745b4457527SToby Isaac 746b4457527SToby Isaac PetscFunctionBegin; 747b4457527SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 748b4457527SToby Isaac PetscValidPointer(intdim, 2); 749b4457527SToby Isaac if (sp->spintdim < 0) { 750b4457527SToby Isaac PetscSection section; 751b4457527SToby Isaac 752b4457527SToby Isaac ierr = PetscDualSpaceGetSection(sp, §ion);CHKERRQ(ierr); 753b4457527SToby Isaac if (section) { 754b4457527SToby Isaac ierr = PetscSectionGetConstrainedStorageSize(section, &(sp->spintdim));CHKERRQ(ierr); 755b4457527SToby Isaac } else sp->spintdim = 0; 756b4457527SToby Isaac } 757b4457527SToby Isaac *intdim = sp->spintdim; 758b4457527SToby Isaac PetscFunctionReturn(0); 759b4457527SToby Isaac } 760b4457527SToby Isaac 761b4457527SToby Isaac /*@ 762b4457527SToby Isaac PetscDualSpaceGetUniform - Whether this dual space is uniform 763b4457527SToby Isaac 764b4457527SToby Isaac Not collective 765b4457527SToby Isaac 766b4457527SToby Isaac Input Parameters: 767b4457527SToby Isaac . sp - A dual space 768b4457527SToby Isaac 769b4457527SToby Isaac Output Parameters: 770b4457527SToby Isaac . uniform - PETSC_TRUE if (a) the dual space is the same for each point in a stratum of the reference DMPlex, and 771b4457527SToby Isaac (b) every symmetry of each point in the reference DMPlex is also a symmetry of the point's dual space. 772b4457527SToby Isaac 773b4457527SToby Isaac Level: advanced 774b4457527SToby Isaac 775b4457527SToby Isaac Note: all of the usual spaces on simplex or tensor-product elements will be uniform, only reference cells 776b4457527SToby Isaac with non-uniform strata (like trianguar-prisms) or anisotropic hp dual spaces will not be uniform. 777b4457527SToby Isaac 778b4457527SToby Isaac .seealso: PetscDualSpaceGetPointSubspace(), PetscDualSpaceGetSymmetries() 779b4457527SToby Isaac @*/ 780b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace sp, PetscBool *uniform) 781b4457527SToby Isaac { 782b4457527SToby Isaac PetscFunctionBegin; 783b4457527SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 784b4457527SToby Isaac PetscValidPointer(uniform, 2); 785b4457527SToby Isaac *uniform = sp->uniform; 786b4457527SToby Isaac PetscFunctionReturn(0); 787b4457527SToby Isaac } 788b4457527SToby Isaac 78920cf1dd8SToby Isaac /*@C 79020cf1dd8SToby Isaac PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension 79120cf1dd8SToby Isaac 79220cf1dd8SToby Isaac Not collective 79320cf1dd8SToby Isaac 79420cf1dd8SToby Isaac Input Parameter: 79520cf1dd8SToby Isaac . sp - The PetscDualSpace 79620cf1dd8SToby Isaac 79720cf1dd8SToby Isaac Output Parameter: 79820cf1dd8SToby Isaac . numDof - An array of length dim+1 which holds the number of dofs for each dimension 79920cf1dd8SToby Isaac 80020cf1dd8SToby Isaac Level: intermediate 80120cf1dd8SToby Isaac 80220cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 80320cf1dd8SToby Isaac @*/ 80420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof) 80520cf1dd8SToby Isaac { 80620cf1dd8SToby Isaac PetscErrorCode ierr; 80720cf1dd8SToby Isaac 80820cf1dd8SToby Isaac PetscFunctionBegin; 80920cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 81020cf1dd8SToby Isaac PetscValidPointer(numDof, 2); 8112c71b3e2SJacob Faibussowitsch PetscCheckFalse(!sp->uniform,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform space does not have a fixed number of dofs for each height"); 812b4457527SToby Isaac if (!sp->numDof) { 813b4457527SToby Isaac DM dm; 814b4457527SToby Isaac PetscInt depth, d; 815b4457527SToby Isaac PetscSection section; 816b4457527SToby Isaac 817b4457527SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 818b4457527SToby Isaac ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 819b4457527SToby Isaac ierr = PetscCalloc1(depth+1,&(sp->numDof));CHKERRQ(ierr); 820b4457527SToby Isaac ierr = PetscDualSpaceGetSection(sp, §ion);CHKERRQ(ierr); 821b4457527SToby Isaac for (d = 0; d <= depth; d++) { 822b4457527SToby Isaac PetscInt dStart, dEnd; 823b4457527SToby Isaac 824b4457527SToby Isaac ierr = DMPlexGetDepthStratum(dm, d, &dStart, &dEnd);CHKERRQ(ierr); 825b4457527SToby Isaac if (dEnd <= dStart) continue; 826b4457527SToby Isaac ierr = PetscSectionGetDof(section, dStart, &(sp->numDof[d]));CHKERRQ(ierr); 827b4457527SToby Isaac 828b4457527SToby Isaac } 829b4457527SToby Isaac } 830b4457527SToby Isaac *numDof = sp->numDof; 8312c71b3e2SJacob Faibussowitsch PetscCheckFalse(!*numDof,PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation"); 83220cf1dd8SToby Isaac PetscFunctionReturn(0); 83320cf1dd8SToby Isaac } 83420cf1dd8SToby Isaac 835b4457527SToby Isaac /* create the section of the right size and set a permutation for topological ordering */ 836b4457527SToby Isaac PetscErrorCode PetscDualSpaceSectionCreate_Internal(PetscDualSpace sp, PetscSection *topSection) 837b4457527SToby Isaac { 838b4457527SToby Isaac DM dm; 839b4457527SToby Isaac PetscInt pStart, pEnd, cStart, cEnd, c, depth, count, i; 840b4457527SToby Isaac PetscInt *seen, *perm; 841b4457527SToby Isaac PetscSection section; 842b4457527SToby Isaac PetscErrorCode ierr; 843b4457527SToby Isaac 844b4457527SToby Isaac PetscFunctionBegin; 845b4457527SToby Isaac dm = sp->dm; 846b4457527SToby Isaac ierr = PetscSectionCreate(PETSC_COMM_SELF, §ion);CHKERRQ(ierr); 847b4457527SToby Isaac ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 848b4457527SToby Isaac ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 849b4457527SToby Isaac ierr = PetscCalloc1(pEnd - pStart, &seen);CHKERRQ(ierr); 850b4457527SToby Isaac ierr = PetscMalloc1(pEnd - pStart, &perm);CHKERRQ(ierr); 851b4457527SToby Isaac ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 852b4457527SToby Isaac ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 853b4457527SToby Isaac for (c = cStart, count = 0; c < cEnd; c++) { 854b4457527SToby Isaac PetscInt closureSize = -1, e; 855b4457527SToby Isaac PetscInt *closure = NULL; 856b4457527SToby Isaac 857b4457527SToby Isaac perm[count++] = c; 858b4457527SToby Isaac seen[c-pStart] = 1; 859b4457527SToby Isaac ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 860b4457527SToby Isaac for (e = 0; e < closureSize; e++) { 861b4457527SToby Isaac PetscInt point = closure[2*e]; 862b4457527SToby Isaac 863b4457527SToby Isaac if (seen[point-pStart]) continue; 864b4457527SToby Isaac perm[count++] = point; 865b4457527SToby Isaac seen[point-pStart] = 1; 866b4457527SToby Isaac } 867b4457527SToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 868b4457527SToby Isaac } 8692c71b3e2SJacob Faibussowitsch PetscCheckFalse(count != pEnd - pStart,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad topological ordering"); 870b4457527SToby Isaac for (i = 0; i < pEnd - pStart; i++) if (perm[i] != i) break; 871b4457527SToby Isaac if (i < pEnd - pStart) { 872b4457527SToby Isaac IS permIS; 873b4457527SToby Isaac 874b4457527SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, perm, PETSC_OWN_POINTER, &permIS);CHKERRQ(ierr); 875b4457527SToby Isaac ierr = ISSetPermutation(permIS);CHKERRQ(ierr); 876b4457527SToby Isaac ierr = PetscSectionSetPermutation(section, permIS);CHKERRQ(ierr); 877b4457527SToby Isaac ierr = ISDestroy(&permIS);CHKERRQ(ierr); 878b4457527SToby Isaac } else { 879b4457527SToby Isaac ierr = PetscFree(perm);CHKERRQ(ierr); 880b4457527SToby Isaac } 881b4457527SToby Isaac ierr = PetscFree(seen);CHKERRQ(ierr); 882b4457527SToby Isaac *topSection = section; 883b4457527SToby Isaac PetscFunctionReturn(0); 884b4457527SToby Isaac } 885b4457527SToby Isaac 886b4457527SToby Isaac /* mark boundary points and set up */ 887b4457527SToby Isaac PetscErrorCode PetscDualSpaceSectionSetUp_Internal(PetscDualSpace sp, PetscSection section) 888b4457527SToby Isaac { 889b4457527SToby Isaac DM dm; 890b4457527SToby Isaac DMLabel boundary; 891b4457527SToby Isaac PetscInt pStart, pEnd, p; 892b4457527SToby Isaac PetscErrorCode ierr; 893b4457527SToby Isaac 894b4457527SToby Isaac PetscFunctionBegin; 895b4457527SToby Isaac dm = sp->dm; 896b4457527SToby Isaac ierr = DMLabelCreate(PETSC_COMM_SELF,"boundary",&boundary);CHKERRQ(ierr); 897b4457527SToby Isaac ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 898b4457527SToby Isaac ierr = DMPlexMarkBoundaryFaces(dm,1,boundary);CHKERRQ(ierr); 899b4457527SToby Isaac ierr = DMPlexLabelComplete(dm,boundary);CHKERRQ(ierr); 900b4457527SToby Isaac ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 901b4457527SToby Isaac for (p = pStart; p < pEnd; p++) { 902b4457527SToby Isaac PetscInt bval; 903b4457527SToby Isaac 904b4457527SToby Isaac ierr = DMLabelGetValue(boundary, p, &bval);CHKERRQ(ierr); 905b4457527SToby Isaac if (bval == 1) { 906b4457527SToby Isaac PetscInt dof; 907b4457527SToby Isaac 908b4457527SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 909b4457527SToby Isaac ierr = PetscSectionSetConstraintDof(section, p, dof);CHKERRQ(ierr); 910b4457527SToby Isaac } 911b4457527SToby Isaac } 912b4457527SToby Isaac ierr = DMLabelDestroy(&boundary);CHKERRQ(ierr); 9131e1ea65dSPierre Jolivet ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 914b4457527SToby Isaac PetscFunctionReturn(0); 915b4457527SToby Isaac } 916b4457527SToby Isaac 917a4ce7ad1SMatthew G. Knepley /*@ 918b4457527SToby Isaac PetscDualSpaceGetSection - Create a PetscSection over the reference cell with the layout from this space 919a4ce7ad1SMatthew G. Knepley 920a4ce7ad1SMatthew G. Knepley Collective on sp 921a4ce7ad1SMatthew G. Knepley 922a4ce7ad1SMatthew G. Knepley Input Parameters: 923f0fc11ceSJed Brown . sp - The PetscDualSpace 924a4ce7ad1SMatthew G. Knepley 925a4ce7ad1SMatthew G. Knepley Output Parameter: 926a4ce7ad1SMatthew G. Knepley . section - The section 927a4ce7ad1SMatthew G. Knepley 928a4ce7ad1SMatthew G. Knepley Level: advanced 929a4ce7ad1SMatthew G. Knepley 930a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate(), DMPLEX 931a4ce7ad1SMatthew G. Knepley @*/ 932b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace sp, PetscSection *section) 93320cf1dd8SToby Isaac { 934b4457527SToby Isaac PetscInt pStart, pEnd, p; 935b4457527SToby Isaac PetscErrorCode ierr; 936b4457527SToby Isaac 937b4457527SToby Isaac PetscFunctionBegin; 938b4457527SToby Isaac if (!sp->pointSection) { 939b4457527SToby Isaac /* mark the boundary */ 940b4457527SToby Isaac ierr = PetscDualSpaceSectionCreate_Internal(sp, &(sp->pointSection));CHKERRQ(ierr); 941b4457527SToby Isaac ierr = DMPlexGetChart(sp->dm,&pStart,&pEnd);CHKERRQ(ierr); 942b4457527SToby Isaac for (p = pStart; p < pEnd; p++) { 943b4457527SToby Isaac PetscDualSpace psp; 944b4457527SToby Isaac 945b4457527SToby Isaac ierr = PetscDualSpaceGetPointSubspace(sp, p, &psp);CHKERRQ(ierr); 946b4457527SToby Isaac if (psp) { 947b4457527SToby Isaac PetscInt dof; 948b4457527SToby Isaac 949b4457527SToby Isaac ierr = PetscDualSpaceGetInteriorDimension(psp, &dof);CHKERRQ(ierr); 950b4457527SToby Isaac ierr = PetscSectionSetDof(sp->pointSection,p,dof);CHKERRQ(ierr); 951b4457527SToby Isaac } 952b4457527SToby Isaac } 953b4457527SToby Isaac ierr = PetscDualSpaceSectionSetUp_Internal(sp,sp->pointSection);CHKERRQ(ierr); 954b4457527SToby Isaac } 955b4457527SToby Isaac *section = sp->pointSection; 956b4457527SToby Isaac PetscFunctionReturn(0); 957b4457527SToby Isaac } 958b4457527SToby Isaac 959b4457527SToby Isaac /* this assumes that all of the point dual spaces store their interior dofs first, which is true when the point DMs 960b4457527SToby Isaac * have one cell */ 961b4457527SToby Isaac PetscErrorCode PetscDualSpacePushForwardSubspaces_Internal(PetscDualSpace sp, PetscInt sStart, PetscInt sEnd) 962b4457527SToby Isaac { 963b4457527SToby Isaac PetscReal *sv0, *v0, *J; 964b4457527SToby Isaac PetscSection section; 965b4457527SToby Isaac PetscInt dim, s, k; 96620cf1dd8SToby Isaac DM dm; 96720cf1dd8SToby Isaac PetscErrorCode ierr; 96820cf1dd8SToby Isaac 96920cf1dd8SToby Isaac PetscFunctionBegin; 97020cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 971b4457527SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 972b4457527SToby Isaac ierr = PetscDualSpaceGetSection(sp, §ion);CHKERRQ(ierr); 973b4457527SToby Isaac ierr = PetscMalloc3(dim, &v0, dim, &sv0, dim*dim, &J);CHKERRQ(ierr); 974b4457527SToby Isaac ierr = PetscDualSpaceGetFormDegree(sp, &k);CHKERRQ(ierr); 975b4457527SToby Isaac for (s = sStart; s < sEnd; s++) { 976b4457527SToby Isaac PetscReal detJ, hdetJ; 977b4457527SToby Isaac PetscDualSpace ssp; 978b4457527SToby Isaac PetscInt dof, off, f, sdim; 979b4457527SToby Isaac PetscInt i, j; 980b4457527SToby Isaac DM sdm; 98120cf1dd8SToby Isaac 982b4457527SToby Isaac ierr = PetscDualSpaceGetPointSubspace(sp, s, &ssp);CHKERRQ(ierr); 983b4457527SToby Isaac if (!ssp) continue; 984b4457527SToby Isaac ierr = PetscSectionGetDof(section, s, &dof);CHKERRQ(ierr); 985b4457527SToby Isaac ierr = PetscSectionGetOffset(section, s, &off);CHKERRQ(ierr); 986b4457527SToby Isaac /* get the first vertex of the reference cell */ 987b4457527SToby Isaac ierr = PetscDualSpaceGetDM(ssp, &sdm);CHKERRQ(ierr); 988b4457527SToby Isaac ierr = DMGetDimension(sdm, &sdim);CHKERRQ(ierr); 989b4457527SToby Isaac ierr = DMPlexComputeCellGeometryAffineFEM(sdm, 0, sv0, NULL, NULL, &hdetJ);CHKERRQ(ierr); 990b4457527SToby Isaac ierr = DMPlexComputeCellGeometryAffineFEM(dm, s, v0, J, NULL, &detJ);CHKERRQ(ierr); 991b4457527SToby Isaac /* compactify Jacobian */ 992b4457527SToby Isaac for (i = 0; i < dim; i++) for (j = 0; j < sdim; j++) J[i* sdim + j] = J[i * dim + j]; 993b4457527SToby Isaac for (f = 0; f < dof; f++) { 994b4457527SToby Isaac PetscQuadrature fn; 99520cf1dd8SToby Isaac 996b4457527SToby Isaac ierr = PetscDualSpaceGetFunctional(ssp, f, &fn);CHKERRQ(ierr); 997b4457527SToby Isaac ierr = PetscQuadraturePushForward(fn, dim, sv0, v0, J, k, &(sp->functional[off+f]));CHKERRQ(ierr); 99820cf1dd8SToby Isaac } 99920cf1dd8SToby Isaac } 1000b4457527SToby Isaac ierr = PetscFree3(v0, sv0, J);CHKERRQ(ierr); 100120cf1dd8SToby Isaac PetscFunctionReturn(0); 100220cf1dd8SToby Isaac } 100320cf1dd8SToby Isaac 100420cf1dd8SToby Isaac /*@C 100520cf1dd8SToby Isaac PetscDualSpaceApply - Apply a functional from the dual space basis to an input function 100620cf1dd8SToby Isaac 100720cf1dd8SToby Isaac Input Parameters: 100820cf1dd8SToby Isaac + sp - The PetscDualSpace object 100920cf1dd8SToby Isaac . f - The basis functional index 101020cf1dd8SToby Isaac . time - The time 101120cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian) (or evaluated at the coordinates of the functional) 101220cf1dd8SToby Isaac . numComp - The number of components for the function 101320cf1dd8SToby Isaac . func - The input function 101420cf1dd8SToby Isaac - ctx - A context for the function 101520cf1dd8SToby Isaac 101620cf1dd8SToby Isaac Output Parameter: 101720cf1dd8SToby Isaac . value - numComp output values 101820cf1dd8SToby Isaac 101920cf1dd8SToby Isaac Note: The calling sequence for the callback func is given by: 102020cf1dd8SToby Isaac 102120cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[], 102220cf1dd8SToby Isaac $ PetscInt numComponents, PetscScalar values[], void *ctx) 102320cf1dd8SToby Isaac 1024a4ce7ad1SMatthew G. Knepley Level: beginner 102520cf1dd8SToby Isaac 102620cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate() 102720cf1dd8SToby Isaac @*/ 102820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApply(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt numComp, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value) 102920cf1dd8SToby Isaac { 103020cf1dd8SToby Isaac PetscErrorCode ierr; 103120cf1dd8SToby Isaac 103220cf1dd8SToby Isaac PetscFunctionBegin; 103320cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 103420cf1dd8SToby Isaac PetscValidPointer(cgeom, 4); 103520cf1dd8SToby Isaac PetscValidPointer(value, 8); 103620cf1dd8SToby Isaac ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr); 103720cf1dd8SToby Isaac PetscFunctionReturn(0); 103820cf1dd8SToby Isaac } 103920cf1dd8SToby Isaac 104020cf1dd8SToby Isaac /*@C 1041b4457527SToby Isaac PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllData() 104220cf1dd8SToby Isaac 104320cf1dd8SToby Isaac Input Parameters: 104420cf1dd8SToby Isaac + sp - The PetscDualSpace object 1045b4457527SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllData() 104620cf1dd8SToby Isaac 104720cf1dd8SToby Isaac Output Parameter: 104820cf1dd8SToby Isaac . spValue - The values of all dual space functionals 104920cf1dd8SToby Isaac 1050a4ce7ad1SMatthew G. Knepley Level: beginner 105120cf1dd8SToby Isaac 105220cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate() 105320cf1dd8SToby Isaac @*/ 105420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 105520cf1dd8SToby Isaac { 105620cf1dd8SToby Isaac PetscErrorCode ierr; 105720cf1dd8SToby Isaac 105820cf1dd8SToby Isaac PetscFunctionBegin; 105920cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 106020cf1dd8SToby Isaac ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr); 106120cf1dd8SToby Isaac PetscFunctionReturn(0); 106220cf1dd8SToby Isaac } 106320cf1dd8SToby Isaac 106420cf1dd8SToby Isaac /*@C 1065b4457527SToby Isaac PetscDualSpaceApplyInterior - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetInteriorData() 1066b4457527SToby Isaac 1067b4457527SToby Isaac Input Parameters: 1068b4457527SToby Isaac + sp - The PetscDualSpace object 1069b4457527SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetInteriorData() 1070b4457527SToby Isaac 1071b4457527SToby Isaac Output Parameter: 1072b4457527SToby Isaac . spValue - The values of interior dual space functionals 1073b4457527SToby Isaac 1074b4457527SToby Isaac Level: beginner 1075b4457527SToby Isaac 1076b4457527SToby Isaac .seealso: PetscDualSpaceCreate() 1077b4457527SToby Isaac @*/ 1078b4457527SToby Isaac PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 1079b4457527SToby Isaac { 1080b4457527SToby Isaac PetscErrorCode ierr; 1081b4457527SToby Isaac 1082b4457527SToby Isaac PetscFunctionBegin; 1083b4457527SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1084b4457527SToby Isaac ierr = (*sp->ops->applyint)(sp, pointEval, spValue);CHKERRQ(ierr); 1085b4457527SToby Isaac PetscFunctionReturn(0); 1086b4457527SToby Isaac } 1087b4457527SToby Isaac 1088b4457527SToby Isaac /*@C 108920cf1dd8SToby Isaac PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional. 109020cf1dd8SToby Isaac 109120cf1dd8SToby Isaac Input Parameters: 109220cf1dd8SToby Isaac + sp - The PetscDualSpace object 109320cf1dd8SToby Isaac . f - The basis functional index 109420cf1dd8SToby Isaac . time - The time 109520cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian) 109620cf1dd8SToby Isaac . Nc - The number of components for the function 109720cf1dd8SToby Isaac . func - The input function 109820cf1dd8SToby Isaac - ctx - A context for the function 109920cf1dd8SToby Isaac 110020cf1dd8SToby Isaac Output Parameter: 110120cf1dd8SToby Isaac . value - The output value 110220cf1dd8SToby Isaac 110320cf1dd8SToby Isaac Note: The calling sequence for the callback func is given by: 110420cf1dd8SToby Isaac 110520cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[], 110620cf1dd8SToby Isaac $ PetscInt numComponents, PetscScalar values[], void *ctx) 110720cf1dd8SToby Isaac 110820cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral 110920cf1dd8SToby Isaac 111020cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x) 111120cf1dd8SToby Isaac 111220cf1dd8SToby Isaac where both n and f have Nc components. 111320cf1dd8SToby Isaac 1114a4ce7ad1SMatthew G. Knepley Level: beginner 111520cf1dd8SToby Isaac 111620cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate() 111720cf1dd8SToby Isaac @*/ 111820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value) 111920cf1dd8SToby Isaac { 112020cf1dd8SToby Isaac DM dm; 112120cf1dd8SToby Isaac PetscQuadrature n; 112220cf1dd8SToby Isaac const PetscReal *points, *weights; 112320cf1dd8SToby Isaac PetscReal x[3]; 112420cf1dd8SToby Isaac PetscScalar *val; 112520cf1dd8SToby Isaac PetscInt dim, dE, qNc, c, Nq, q; 112620cf1dd8SToby Isaac PetscBool isAffine; 112720cf1dd8SToby Isaac PetscErrorCode ierr; 112820cf1dd8SToby Isaac 112920cf1dd8SToby Isaac PetscFunctionBegin; 113020cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1131064a246eSJacob Faibussowitsch PetscValidPointer(value, 8); 113220cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 113320cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 113420cf1dd8SToby Isaac ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 11352c71b3e2SJacob Faibussowitsch PetscCheckFalse(dim != cgeom->dim,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim); 11362c71b3e2SJacob Faibussowitsch PetscCheckFalse(qNc != Nc,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 113720cf1dd8SToby Isaac ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 113820cf1dd8SToby Isaac *value = 0.0; 113920cf1dd8SToby Isaac isAffine = cgeom->isAffine; 114020cf1dd8SToby Isaac dE = cgeom->dimEmbed; 114120cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 114220cf1dd8SToby Isaac if (isAffine) { 114320cf1dd8SToby Isaac CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x); 114420cf1dd8SToby Isaac ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr); 114520cf1dd8SToby Isaac } else { 114620cf1dd8SToby Isaac ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr); 114720cf1dd8SToby Isaac } 114820cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 114920cf1dd8SToby Isaac *value += val[c]*weights[q*Nc+c]; 115020cf1dd8SToby Isaac } 115120cf1dd8SToby Isaac } 115220cf1dd8SToby Isaac ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 115320cf1dd8SToby Isaac PetscFunctionReturn(0); 115420cf1dd8SToby Isaac } 115520cf1dd8SToby Isaac 115620cf1dd8SToby Isaac /*@C 1157b4457527SToby Isaac PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllData() 115820cf1dd8SToby Isaac 115920cf1dd8SToby Isaac Input Parameters: 116020cf1dd8SToby Isaac + sp - The PetscDualSpace object 1161b4457527SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllData() 116220cf1dd8SToby Isaac 116320cf1dd8SToby Isaac Output Parameter: 116420cf1dd8SToby Isaac . spValue - The values of all dual space functionals 116520cf1dd8SToby Isaac 1166a4ce7ad1SMatthew G. Knepley Level: beginner 116720cf1dd8SToby Isaac 116820cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate() 116920cf1dd8SToby Isaac @*/ 117020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 117120cf1dd8SToby Isaac { 1172b4457527SToby Isaac Vec pointValues, dofValues; 1173b4457527SToby Isaac Mat allMat; 117420cf1dd8SToby Isaac PetscErrorCode ierr; 117520cf1dd8SToby Isaac 117620cf1dd8SToby Isaac PetscFunctionBegin; 117720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 117820cf1dd8SToby Isaac PetscValidScalarPointer(pointEval, 2); 1179064a246eSJacob Faibussowitsch PetscValidScalarPointer(spValue, 3); 1180b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp, NULL, &allMat);CHKERRQ(ierr); 1181b4457527SToby Isaac if (!(sp->allNodeValues)) { 1182b4457527SToby Isaac ierr = MatCreateVecs(allMat, &(sp->allNodeValues), NULL);CHKERRQ(ierr); 118320cf1dd8SToby Isaac } 1184b4457527SToby Isaac pointValues = sp->allNodeValues; 1185b4457527SToby Isaac if (!(sp->allDofValues)) { 1186b4457527SToby Isaac ierr = MatCreateVecs(allMat, NULL, &(sp->allDofValues));CHKERRQ(ierr); 118720cf1dd8SToby Isaac } 1188b4457527SToby Isaac dofValues = sp->allDofValues; 1189b4457527SToby Isaac ierr = VecPlaceArray(pointValues, pointEval);CHKERRQ(ierr); 1190b4457527SToby Isaac ierr = VecPlaceArray(dofValues, spValue);CHKERRQ(ierr); 1191b4457527SToby Isaac ierr = MatMult(allMat, pointValues, dofValues);CHKERRQ(ierr); 1192b4457527SToby Isaac ierr = VecResetArray(dofValues);CHKERRQ(ierr); 1193b4457527SToby Isaac ierr = VecResetArray(pointValues);CHKERRQ(ierr); 1194b4457527SToby Isaac PetscFunctionReturn(0); 119520cf1dd8SToby Isaac } 1196b4457527SToby Isaac 1197b4457527SToby Isaac /*@C 1198b4457527SToby Isaac PetscDualSpaceApplyInteriorDefault - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetInteriorData() 1199b4457527SToby Isaac 1200b4457527SToby Isaac Input Parameters: 1201b4457527SToby Isaac + sp - The PetscDualSpace object 1202b4457527SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetInteriorData() 1203b4457527SToby Isaac 1204b4457527SToby Isaac Output Parameter: 1205b4457527SToby Isaac . spValue - The values of interior dual space functionals 1206b4457527SToby Isaac 1207b4457527SToby Isaac Level: beginner 1208b4457527SToby Isaac 1209b4457527SToby Isaac .seealso: PetscDualSpaceCreate() 1210b4457527SToby Isaac @*/ 1211b4457527SToby Isaac PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 1212b4457527SToby Isaac { 1213b4457527SToby Isaac Vec pointValues, dofValues; 1214b4457527SToby Isaac Mat intMat; 1215b4457527SToby Isaac PetscErrorCode ierr; 1216b4457527SToby Isaac 1217b4457527SToby Isaac PetscFunctionBegin; 1218b4457527SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1219b4457527SToby Isaac PetscValidScalarPointer(pointEval, 2); 1220064a246eSJacob Faibussowitsch PetscValidScalarPointer(spValue, 3); 1221b4457527SToby Isaac ierr = PetscDualSpaceGetInteriorData(sp, NULL, &intMat);CHKERRQ(ierr); 1222b4457527SToby Isaac if (!(sp->intNodeValues)) { 1223b4457527SToby Isaac ierr = MatCreateVecs(intMat, &(sp->intNodeValues), NULL);CHKERRQ(ierr); 1224b4457527SToby Isaac } 1225b4457527SToby Isaac pointValues = sp->intNodeValues; 1226b4457527SToby Isaac if (!(sp->intDofValues)) { 1227b4457527SToby Isaac ierr = MatCreateVecs(intMat, NULL, &(sp->intDofValues));CHKERRQ(ierr); 1228b4457527SToby Isaac } 1229b4457527SToby Isaac dofValues = sp->intDofValues; 1230b4457527SToby Isaac ierr = VecPlaceArray(pointValues, pointEval);CHKERRQ(ierr); 1231b4457527SToby Isaac ierr = VecPlaceArray(dofValues, spValue);CHKERRQ(ierr); 1232b4457527SToby Isaac ierr = MatMult(intMat, pointValues, dofValues);CHKERRQ(ierr); 1233b4457527SToby Isaac ierr = VecResetArray(dofValues);CHKERRQ(ierr); 1234b4457527SToby Isaac ierr = VecResetArray(pointValues);CHKERRQ(ierr); 123520cf1dd8SToby Isaac PetscFunctionReturn(0); 123620cf1dd8SToby Isaac } 123720cf1dd8SToby Isaac 1238a4ce7ad1SMatthew G. Knepley /*@ 1239b4457527SToby Isaac PetscDualSpaceGetAllData - Get all quadrature nodes from this space, and the matrix that sends quadrature node values to degree-of-freedom values 1240a4ce7ad1SMatthew G. Knepley 1241a4ce7ad1SMatthew G. Knepley Input Parameter: 1242a4ce7ad1SMatthew G. Knepley . sp - The dualspace 1243a4ce7ad1SMatthew G. Knepley 1244d8d19677SJose E. Roman Output Parameters: 1245b4457527SToby Isaac + allNodes - A PetscQuadrature object containing all evaluation nodes 1246b4457527SToby Isaac - allMat - A Mat for the node-to-dof transformation 1247a4ce7ad1SMatthew G. Knepley 1248a4ce7ad1SMatthew G. Knepley Level: advanced 1249a4ce7ad1SMatthew G. Knepley 1250a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate() 1251a4ce7ad1SMatthew G. Knepley @*/ 1252b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat) 125320cf1dd8SToby Isaac { 125420cf1dd8SToby Isaac PetscErrorCode ierr; 125520cf1dd8SToby Isaac 125620cf1dd8SToby Isaac PetscFunctionBegin; 125720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1258b4457527SToby Isaac if (allNodes) PetscValidPointer(allNodes,2); 1259b4457527SToby Isaac if (allMat) PetscValidPointer(allMat,3); 1260b4457527SToby Isaac if ((!sp->allNodes || !sp->allMat) && sp->ops->createalldata) { 1261b4457527SToby Isaac PetscQuadrature qpoints; 1262b4457527SToby Isaac Mat amat; 1263b4457527SToby Isaac 1264b4457527SToby Isaac ierr = (*sp->ops->createalldata)(sp,&qpoints,&amat);CHKERRQ(ierr); 1265b4457527SToby Isaac ierr = PetscQuadratureDestroy(&(sp->allNodes));CHKERRQ(ierr); 1266b4457527SToby Isaac ierr = MatDestroy(&(sp->allMat));CHKERRQ(ierr); 1267b4457527SToby Isaac sp->allNodes = qpoints; 1268b4457527SToby Isaac sp->allMat = amat; 126920cf1dd8SToby Isaac } 1270b4457527SToby Isaac if (allNodes) *allNodes = sp->allNodes; 1271b4457527SToby Isaac if (allMat) *allMat = sp->allMat; 127220cf1dd8SToby Isaac PetscFunctionReturn(0); 127320cf1dd8SToby Isaac } 127420cf1dd8SToby Isaac 1275a4ce7ad1SMatthew G. Knepley /*@ 1276b4457527SToby Isaac PetscDualSpaceCreateAllDataDefault - Create all evaluation nodes and the node-to-dof matrix by examining functionals 1277a4ce7ad1SMatthew G. Knepley 1278a4ce7ad1SMatthew G. Knepley Input Parameter: 1279a4ce7ad1SMatthew G. Knepley . sp - The dualspace 1280a4ce7ad1SMatthew G. Knepley 1281d8d19677SJose E. Roman Output Parameters: 1282b4457527SToby Isaac + allNodes - A PetscQuadrature object containing all evaluation nodes 1283b4457527SToby Isaac - allMat - A Mat for the node-to-dof transformation 1284a4ce7ad1SMatthew G. Knepley 1285a4ce7ad1SMatthew G. Knepley Level: advanced 1286a4ce7ad1SMatthew G. Knepley 1287a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate() 1288a4ce7ad1SMatthew G. Knepley @*/ 1289b4457527SToby Isaac PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat) 129020cf1dd8SToby Isaac { 129120cf1dd8SToby Isaac PetscInt spdim; 129220cf1dd8SToby Isaac PetscInt numPoints, offset; 129320cf1dd8SToby Isaac PetscReal *points; 129420cf1dd8SToby Isaac PetscInt f, dim; 1295b4457527SToby Isaac PetscInt Nc, nrows, ncols; 1296b4457527SToby Isaac PetscInt maxNumPoints; 129720cf1dd8SToby Isaac PetscQuadrature q; 1298b4457527SToby Isaac Mat A; 129920cf1dd8SToby Isaac PetscErrorCode ierr; 130020cf1dd8SToby Isaac 130120cf1dd8SToby Isaac PetscFunctionBegin; 1302b4457527SToby Isaac ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 130320cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr); 130420cf1dd8SToby Isaac if (!spdim) { 1305b4457527SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allNodes);CHKERRQ(ierr); 1306b4457527SToby Isaac ierr = PetscQuadratureSetData(*allNodes,0,0,0,NULL,NULL);CHKERRQ(ierr); 130720cf1dd8SToby Isaac } 1308b4457527SToby Isaac nrows = spdim; 130920cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr); 131020cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr); 1311b4457527SToby Isaac maxNumPoints = numPoints; 131220cf1dd8SToby Isaac for (f = 1; f < spdim; f++) { 131320cf1dd8SToby Isaac PetscInt Np; 131420cf1dd8SToby Isaac 131520cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr); 131620cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr); 131720cf1dd8SToby Isaac numPoints += Np; 1318b4457527SToby Isaac maxNumPoints = PetscMax(maxNumPoints,Np); 131920cf1dd8SToby Isaac } 1320b4457527SToby Isaac ncols = numPoints * Nc; 132120cf1dd8SToby Isaac ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 1322b4457527SToby Isaac ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, nrows, ncols, maxNumPoints * Nc, NULL, &A);CHKERRQ(ierr); 132320cf1dd8SToby Isaac for (f = 0, offset = 0; f < spdim; f++) { 1324b4457527SToby Isaac const PetscReal *p, *w; 132520cf1dd8SToby Isaac PetscInt Np, i; 1326b4457527SToby Isaac PetscInt fnc; 132720cf1dd8SToby Isaac 132820cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr); 1329b4457527SToby Isaac ierr = PetscQuadratureGetData(q,NULL,&fnc,&Np,&p,&w);CHKERRQ(ierr); 13302c71b3e2SJacob Faibussowitsch PetscCheckFalse(fnc != Nc,PETSC_COMM_SELF, PETSC_ERR_PLIB, "functional component mismatch"); 1331b4457527SToby Isaac for (i = 0; i < Np * dim; i++) { 1332b4457527SToby Isaac points[offset* dim + i] = p[i]; 1333b4457527SToby Isaac } 1334b4457527SToby Isaac for (i = 0; i < Np * Nc; i++) { 1335b4457527SToby Isaac ierr = MatSetValue(A, f, offset * Nc, w[i], INSERT_VALUES);CHKERRQ(ierr); 1336b4457527SToby Isaac } 1337b4457527SToby Isaac offset += Np; 1338b4457527SToby Isaac } 1339b4457527SToby Isaac ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1340b4457527SToby Isaac ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1341b4457527SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allNodes);CHKERRQ(ierr); 1342b4457527SToby Isaac ierr = PetscQuadratureSetData(*allNodes,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 1343b4457527SToby Isaac *allMat = A; 1344b4457527SToby Isaac PetscFunctionReturn(0); 1345b4457527SToby Isaac } 1346b4457527SToby Isaac 1347b4457527SToby Isaac /*@ 1348b4457527SToby Isaac PetscDualSpaceGetInteriorData - Get all quadrature points necessary to compute the interior degrees of freedom from 1349b4457527SToby Isaac this space, as well as the matrix that computes the degrees of freedom from the quadrature values. Degrees of 1350b4457527SToby Isaac freedom are interior degrees of freedom if they belong (by PetscDualSpaceGetSection()) to interior points in the 1351b4457527SToby Isaac reference DMPlex: complementary boundary degrees of freedom are marked as constrained in the section returned by 1352b4457527SToby Isaac PetscDualSpaceGetSection()). 1353b4457527SToby Isaac 1354b4457527SToby Isaac Input Parameter: 1355b4457527SToby Isaac . sp - The dualspace 1356b4457527SToby Isaac 1357d8d19677SJose E. Roman Output Parameters: 1358b4457527SToby Isaac + intNodes - A PetscQuadrature object containing all evaluation points needed to evaluate interior degrees of freedom 1359b4457527SToby Isaac - intMat - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is 1360b4457527SToby Isaac the size of the constrained layout (PetscSectionGetConstrainStorageSize()) of the dual space section, 1361b4457527SToby Isaac npoints is the number of points in intNodes and nc is PetscDualSpaceGetNumComponents(). 1362b4457527SToby Isaac 1363b4457527SToby Isaac Level: advanced 1364b4457527SToby Isaac 1365b4457527SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceGetDimension(), PetscDualSpaceGetNumComponents(), PetscQuadratureGetData() 1366b4457527SToby Isaac @*/ 1367b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat) 1368b4457527SToby Isaac { 1369b4457527SToby Isaac PetscErrorCode ierr; 1370b4457527SToby Isaac 1371b4457527SToby Isaac PetscFunctionBegin; 1372b4457527SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1373b4457527SToby Isaac if (intNodes) PetscValidPointer(intNodes,2); 1374b4457527SToby Isaac if (intMat) PetscValidPointer(intMat,3); 1375b4457527SToby Isaac if ((!sp->intNodes || !sp->intMat) && sp->ops->createintdata) { 1376b4457527SToby Isaac PetscQuadrature qpoints; 1377b4457527SToby Isaac Mat imat; 1378b4457527SToby Isaac 1379b4457527SToby Isaac ierr = (*sp->ops->createintdata)(sp,&qpoints,&imat);CHKERRQ(ierr); 1380b4457527SToby Isaac ierr = PetscQuadratureDestroy(&(sp->intNodes));CHKERRQ(ierr); 1381b4457527SToby Isaac ierr = MatDestroy(&(sp->intMat));CHKERRQ(ierr); 1382b4457527SToby Isaac sp->intNodes = qpoints; 1383b4457527SToby Isaac sp->intMat = imat; 1384b4457527SToby Isaac } 1385b4457527SToby Isaac if (intNodes) *intNodes = sp->intNodes; 1386b4457527SToby Isaac if (intMat) *intMat = sp->intMat; 1387b4457527SToby Isaac PetscFunctionReturn(0); 1388b4457527SToby Isaac } 1389b4457527SToby Isaac 1390b4457527SToby Isaac /*@ 1391b4457527SToby Isaac PetscDualSpaceCreateInteriorDataDefault - Create quadrature points by examining interior functionals and create the matrix mapping quadrature point values to interior dual space values 1392b4457527SToby Isaac 1393b4457527SToby Isaac Input Parameter: 1394b4457527SToby Isaac . sp - The dualspace 1395b4457527SToby Isaac 1396d8d19677SJose E. Roman Output Parameters: 1397b4457527SToby Isaac + intNodes - A PetscQuadrature object containing all evaluation points needed to evaluate interior degrees of freedom 1398b4457527SToby Isaac - intMat - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is 1399b4457527SToby Isaac the size of the constrained layout (PetscSectionGetConstrainStorageSize()) of the dual space section, 1400b4457527SToby Isaac npoints is the number of points in allNodes and nc is PetscDualSpaceGetNumComponents(). 1401b4457527SToby Isaac 1402b4457527SToby Isaac Level: advanced 1403b4457527SToby Isaac 1404b4457527SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceGetInteriorData() 1405b4457527SToby Isaac @*/ 1406b4457527SToby Isaac PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat) 1407b4457527SToby Isaac { 1408b4457527SToby Isaac DM dm; 1409b4457527SToby Isaac PetscInt spdim0; 1410b4457527SToby Isaac PetscInt Nc; 1411b4457527SToby Isaac PetscInt pStart, pEnd, p, f; 1412b4457527SToby Isaac PetscSection section; 1413b4457527SToby Isaac PetscInt numPoints, offset, matoffset; 1414b4457527SToby Isaac PetscReal *points; 1415b4457527SToby Isaac PetscInt dim; 1416b4457527SToby Isaac PetscInt *nnz; 1417b4457527SToby Isaac PetscQuadrature q; 1418b4457527SToby Isaac Mat imat; 1419b4457527SToby Isaac PetscErrorCode ierr; 1420b4457527SToby Isaac 1421b4457527SToby Isaac PetscFunctionBegin; 1422b4457527SToby Isaac PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1); 1423b4457527SToby Isaac ierr = PetscDualSpaceGetSection(sp, §ion);CHKERRQ(ierr); 1424b4457527SToby Isaac ierr = PetscSectionGetConstrainedStorageSize(section, &spdim0);CHKERRQ(ierr); 1425b4457527SToby Isaac if (!spdim0) { 1426b4457527SToby Isaac *intNodes = NULL; 1427b4457527SToby Isaac *intMat = NULL; 1428b4457527SToby Isaac PetscFunctionReturn(0); 1429b4457527SToby Isaac } 1430b4457527SToby Isaac ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 1431b4457527SToby Isaac ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 1432b4457527SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 1433b4457527SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1434b4457527SToby Isaac ierr = PetscMalloc1(spdim0, &nnz);CHKERRQ(ierr); 1435b4457527SToby Isaac for (p = pStart, f = 0, numPoints = 0; p < pEnd; p++) { 1436b4457527SToby Isaac PetscInt dof, cdof, off, d; 1437b4457527SToby Isaac 1438b4457527SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 1439b4457527SToby Isaac ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 1440b4457527SToby Isaac if (!(dof - cdof)) continue; 1441b4457527SToby Isaac ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 1442b4457527SToby Isaac for (d = 0; d < dof; d++, off++, f++) { 1443b4457527SToby Isaac PetscInt Np; 1444b4457527SToby Isaac 1445b4457527SToby Isaac ierr = PetscDualSpaceGetFunctional(sp,off,&q);CHKERRQ(ierr); 1446b4457527SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr); 1447b4457527SToby Isaac nnz[f] = Np * Nc; 1448b4457527SToby Isaac numPoints += Np; 1449b4457527SToby Isaac } 1450b4457527SToby Isaac } 1451b4457527SToby Isaac ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, spdim0, numPoints * Nc, 0, nnz, &imat);CHKERRQ(ierr); 1452b4457527SToby Isaac ierr = PetscFree(nnz);CHKERRQ(ierr); 1453b4457527SToby Isaac ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 1454b4457527SToby Isaac for (p = pStart, f = 0, offset = 0, matoffset = 0; p < pEnd; p++) { 1455b4457527SToby Isaac PetscInt dof, cdof, off, d; 1456b4457527SToby Isaac 1457b4457527SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 1458b4457527SToby Isaac ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 1459b4457527SToby Isaac if (!(dof - cdof)) continue; 1460b4457527SToby Isaac ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 1461b4457527SToby Isaac for (d = 0; d < dof; d++, off++, f++) { 1462b4457527SToby Isaac const PetscReal *p; 1463b4457527SToby Isaac const PetscReal *w; 1464b4457527SToby Isaac PetscInt Np, i; 1465b4457527SToby Isaac 1466b4457527SToby Isaac ierr = PetscDualSpaceGetFunctional(sp,off,&q);CHKERRQ(ierr); 1467b4457527SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,&w);CHKERRQ(ierr); 146820cf1dd8SToby Isaac for (i = 0; i < Np * dim; i++) { 146920cf1dd8SToby Isaac points[offset + i] = p[i]; 147020cf1dd8SToby Isaac } 1471b4457527SToby Isaac for (i = 0; i < Np * Nc; i++) { 1472b4457527SToby Isaac ierr = MatSetValue(imat, f, matoffset + i, w[i],INSERT_VALUES);CHKERRQ(ierr); 147320cf1dd8SToby Isaac } 1474b4457527SToby Isaac offset += Np * dim; 1475b4457527SToby Isaac matoffset += Np * Nc; 1476b4457527SToby Isaac } 1477b4457527SToby Isaac } 1478b4457527SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,intNodes);CHKERRQ(ierr); 1479b4457527SToby Isaac ierr = PetscQuadratureSetData(*intNodes,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 1480b4457527SToby Isaac ierr = MatAssemblyBegin(imat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1481b4457527SToby Isaac ierr = MatAssemblyEnd(imat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1482b4457527SToby Isaac *intMat = imat; 148320cf1dd8SToby Isaac PetscFunctionReturn(0); 148420cf1dd8SToby Isaac } 148520cf1dd8SToby Isaac 1486*4f9ab2b4SJed Brown /*@ 1487*4f9ab2b4SJed Brown PetscDualSpaceEqual - Determine if a dual space is equivalent 1488*4f9ab2b4SJed Brown 1489*4f9ab2b4SJed Brown Input Parameters: 1490*4f9ab2b4SJed Brown + A - A PetscDualSpace object 1491*4f9ab2b4SJed Brown - B - Another PetscDualSpace object 1492*4f9ab2b4SJed Brown 1493*4f9ab2b4SJed Brown Output Parameter: 1494*4f9ab2b4SJed Brown . equal - PETSC_TRUE if the dual spaces are equivalent 1495*4f9ab2b4SJed Brown 1496*4f9ab2b4SJed Brown Level: advanced 1497*4f9ab2b4SJed Brown 1498*4f9ab2b4SJed Brown .seealso: PetscDualSpaceCreate() 1499*4f9ab2b4SJed Brown @*/ 1500*4f9ab2b4SJed Brown PetscErrorCode PetscDualSpaceEqual(PetscDualSpace A, PetscDualSpace B, PetscBool *equal) 1501*4f9ab2b4SJed Brown { 1502*4f9ab2b4SJed Brown PetscErrorCode ierr; 1503*4f9ab2b4SJed Brown PetscInt sizeA, sizeB, dimA, dimB; 1504*4f9ab2b4SJed Brown const PetscInt *dofA, *dofB; 1505*4f9ab2b4SJed Brown PetscQuadrature quadA, quadB; 1506*4f9ab2b4SJed Brown Mat matA, matB; 1507*4f9ab2b4SJed Brown 1508*4f9ab2b4SJed Brown PetscFunctionBegin; 1509*4f9ab2b4SJed Brown PetscValidHeaderSpecific(A,PETSCDUALSPACE_CLASSID,1); 1510*4f9ab2b4SJed Brown PetscValidHeaderSpecific(B,PETSCDUALSPACE_CLASSID,2); 1511*4f9ab2b4SJed Brown PetscValidBoolPointer(equal,3); 1512*4f9ab2b4SJed Brown *equal = PETSC_FALSE; 1513*4f9ab2b4SJed Brown ierr = PetscDualSpaceGetDimension(A, &sizeA);CHKERRQ(ierr); 1514*4f9ab2b4SJed Brown ierr = PetscDualSpaceGetDimension(B, &sizeB);CHKERRQ(ierr); 1515*4f9ab2b4SJed Brown if (sizeB != sizeA) { 1516*4f9ab2b4SJed Brown PetscFunctionReturn(0); 1517*4f9ab2b4SJed Brown } 1518*4f9ab2b4SJed Brown ierr = DMGetDimension(A->dm, &dimA);CHKERRQ(ierr); 1519*4f9ab2b4SJed Brown ierr = DMGetDimension(B->dm, &dimB);CHKERRQ(ierr); 1520*4f9ab2b4SJed Brown if (dimA != dimB) { 1521*4f9ab2b4SJed Brown PetscFunctionReturn(0); 1522*4f9ab2b4SJed Brown } 1523*4f9ab2b4SJed Brown 1524*4f9ab2b4SJed Brown ierr = PetscDualSpaceGetNumDof(A, &dofA);CHKERRQ(ierr); 1525*4f9ab2b4SJed Brown ierr = PetscDualSpaceGetNumDof(B, &dofB);CHKERRQ(ierr); 1526*4f9ab2b4SJed Brown for (PetscInt d=0; d<dimA; d++) { 1527*4f9ab2b4SJed Brown if (dofA[d] != dofB[d]) { 1528*4f9ab2b4SJed Brown PetscFunctionReturn(0); 1529*4f9ab2b4SJed Brown } 1530*4f9ab2b4SJed Brown } 1531*4f9ab2b4SJed Brown 1532*4f9ab2b4SJed Brown ierr = PetscDualSpaceGetInteriorData(A, &quadA, &matA);CHKERRQ(ierr); 1533*4f9ab2b4SJed Brown ierr = PetscDualSpaceGetInteriorData(B, &quadB, &matB);CHKERRQ(ierr); 1534*4f9ab2b4SJed Brown if (!quadA && !quadB) { 1535*4f9ab2b4SJed Brown *equal = PETSC_TRUE; 1536*4f9ab2b4SJed Brown } else if (quadA && quadB) { 1537*4f9ab2b4SJed Brown ierr = PetscQuadratureEqual(quadA, quadB, equal);CHKERRQ(ierr); 1538*4f9ab2b4SJed Brown if (*equal == PETSC_FALSE) PetscFunctionReturn(0); 1539*4f9ab2b4SJed Brown if (!matA && !matB) PetscFunctionReturn(0); 1540*4f9ab2b4SJed Brown if (matA && matB) {ierr = MatEqual(matA, matB, equal);CHKERRQ(ierr);} 1541*4f9ab2b4SJed Brown else *equal = PETSC_FALSE; 1542*4f9ab2b4SJed Brown } 1543*4f9ab2b4SJed Brown PetscFunctionReturn(0); 1544*4f9ab2b4SJed Brown } 1545*4f9ab2b4SJed Brown 154620cf1dd8SToby Isaac /*@C 154720cf1dd8SToby Isaac PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid. 154820cf1dd8SToby Isaac 154920cf1dd8SToby Isaac Input Parameters: 155020cf1dd8SToby Isaac + sp - The PetscDualSpace object 155120cf1dd8SToby Isaac . f - The basis functional index 155220cf1dd8SToby Isaac . time - The time 155320cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we currently just use the centroid 155420cf1dd8SToby Isaac . Nc - The number of components for the function 155520cf1dd8SToby Isaac . func - The input function 155620cf1dd8SToby Isaac - ctx - A context for the function 155720cf1dd8SToby Isaac 155820cf1dd8SToby Isaac Output Parameter: 155920cf1dd8SToby Isaac . value - The output value (scalar) 156020cf1dd8SToby Isaac 156120cf1dd8SToby Isaac Note: The calling sequence for the callback func is given by: 156220cf1dd8SToby Isaac 156320cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[], 156420cf1dd8SToby Isaac $ PetscInt numComponents, PetscScalar values[], void *ctx) 156520cf1dd8SToby Isaac 156620cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral 156720cf1dd8SToby Isaac 156820cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x) 156920cf1dd8SToby Isaac 157020cf1dd8SToby Isaac where both n and f have Nc components. 157120cf1dd8SToby Isaac 1572a4ce7ad1SMatthew G. Knepley Level: beginner 157320cf1dd8SToby Isaac 157420cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate() 157520cf1dd8SToby Isaac @*/ 157620cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyFVM(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFVCellGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value) 157720cf1dd8SToby Isaac { 157820cf1dd8SToby Isaac DM dm; 157920cf1dd8SToby Isaac PetscQuadrature n; 158020cf1dd8SToby Isaac const PetscReal *points, *weights; 158120cf1dd8SToby Isaac PetscScalar *val; 158220cf1dd8SToby Isaac PetscInt dimEmbed, qNc, c, Nq, q; 158320cf1dd8SToby Isaac PetscErrorCode ierr; 158420cf1dd8SToby Isaac 158520cf1dd8SToby Isaac PetscFunctionBegin; 158620cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1587064a246eSJacob Faibussowitsch PetscValidPointer(value, 8); 158820cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 158920cf1dd8SToby Isaac ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 159020cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 159120cf1dd8SToby Isaac ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 15922c71b3e2SJacob Faibussowitsch PetscCheckFalse(qNc != Nc,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 159320cf1dd8SToby Isaac ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 159420cf1dd8SToby Isaac *value = 0.; 159520cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 159620cf1dd8SToby Isaac ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr); 159720cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 159820cf1dd8SToby Isaac *value += val[c]*weights[q*Nc+c]; 159920cf1dd8SToby Isaac } 160020cf1dd8SToby Isaac } 160120cf1dd8SToby Isaac ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 160220cf1dd8SToby Isaac PetscFunctionReturn(0); 160320cf1dd8SToby Isaac } 160420cf1dd8SToby Isaac 160520cf1dd8SToby Isaac /*@ 160620cf1dd8SToby Isaac PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a 160720cf1dd8SToby Isaac given height. This assumes that the reference cell is symmetric over points of this height. 160820cf1dd8SToby Isaac 160920cf1dd8SToby Isaac If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and 161020cf1dd8SToby Isaac pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not 161120cf1dd8SToby Isaac support extracting subspaces, then NULL is returned. 161220cf1dd8SToby Isaac 161320cf1dd8SToby Isaac This does not increment the reference count on the returned dual space, and the user should not destroy it. 161420cf1dd8SToby Isaac 161520cf1dd8SToby Isaac Not collective 161620cf1dd8SToby Isaac 161720cf1dd8SToby Isaac Input Parameters: 161820cf1dd8SToby Isaac + sp - the PetscDualSpace object 161920cf1dd8SToby Isaac - height - the height of the mesh point for which the subspace is desired 162020cf1dd8SToby Isaac 162120cf1dd8SToby Isaac Output Parameter: 162220cf1dd8SToby Isaac . subsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the 162320cf1dd8SToby Isaac point, which will be of lesser dimension if height > 0. 162420cf1dd8SToby Isaac 162520cf1dd8SToby Isaac Level: advanced 162620cf1dd8SToby Isaac 162720cf1dd8SToby Isaac .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace 162820cf1dd8SToby Isaac @*/ 162920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp) 163020cf1dd8SToby Isaac { 1631b4457527SToby Isaac PetscInt depth = -1, cStart, cEnd; 1632b4457527SToby Isaac DM dm; 163320cf1dd8SToby Isaac PetscErrorCode ierr; 163420cf1dd8SToby Isaac 163520cf1dd8SToby Isaac PetscFunctionBegin; 163620cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1637064a246eSJacob Faibussowitsch PetscValidPointer(subsp,3); 16382c71b3e2SJacob Faibussowitsch PetscCheckFalse(!(sp->uniform),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform dual space does not have a single dual space at each height"); 163920cf1dd8SToby Isaac *subsp = NULL; 1640b4457527SToby Isaac dm = sp->dm; 1641b4457527SToby Isaac ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 16422c71b3e2SJacob Faibussowitsch PetscCheckFalse(height < 0 || height > depth,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid height"); 1643b4457527SToby Isaac ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 1644b4457527SToby Isaac if (height == 0 && cEnd == cStart + 1) { 1645b4457527SToby Isaac *subsp = sp; 1646b4457527SToby Isaac PetscFunctionReturn(0); 1647b4457527SToby Isaac } 1648b4457527SToby Isaac if (!sp->heightSpaces) { 1649b4457527SToby Isaac PetscInt h; 1650b4457527SToby Isaac ierr = PetscCalloc1(depth+1, &(sp->heightSpaces));CHKERRQ(ierr); 1651b4457527SToby Isaac 1652b4457527SToby Isaac for (h = 0; h <= depth; h++) { 1653b4457527SToby Isaac if (h == 0 && cEnd == cStart + 1) continue; 1654b4457527SToby Isaac if (sp->ops->createheightsubspace) {ierr = (*sp->ops->createheightsubspace)(sp,height,&(sp->heightSpaces[h]));CHKERRQ(ierr);} 1655b4457527SToby Isaac else if (sp->pointSpaces) { 1656b4457527SToby Isaac PetscInt hStart, hEnd; 1657b4457527SToby Isaac 1658b4457527SToby Isaac ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr); 1659b4457527SToby Isaac if (hEnd > hStart) { 1660665f567fSMatthew G. Knepley const char *name; 1661665f567fSMatthew G. Knepley 1662b4457527SToby Isaac ierr = PetscObjectReference((PetscObject)(sp->pointSpaces[hStart]));CHKERRQ(ierr); 1663665f567fSMatthew G. Knepley if (sp->pointSpaces[hStart]) { 1664665f567fSMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) sp, &name);CHKERRQ(ierr); 1665665f567fSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) sp->pointSpaces[hStart], name);CHKERRQ(ierr); 1666665f567fSMatthew G. Knepley } 1667b4457527SToby Isaac sp->heightSpaces[h] = sp->pointSpaces[hStart]; 1668b4457527SToby Isaac } 1669b4457527SToby Isaac } 1670b4457527SToby Isaac } 1671b4457527SToby Isaac } 1672b4457527SToby Isaac *subsp = sp->heightSpaces[height]; 167320cf1dd8SToby Isaac PetscFunctionReturn(0); 167420cf1dd8SToby Isaac } 167520cf1dd8SToby Isaac 167620cf1dd8SToby Isaac /*@ 167720cf1dd8SToby Isaac PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point. 167820cf1dd8SToby Isaac 167920cf1dd8SToby Isaac If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not 168020cf1dd8SToby Isaac defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting 168120cf1dd8SToby Isaac subspaces, then NULL is returned. 168220cf1dd8SToby Isaac 168320cf1dd8SToby Isaac This does not increment the reference count on the returned dual space, and the user should not destroy it. 168420cf1dd8SToby Isaac 168520cf1dd8SToby Isaac Not collective 168620cf1dd8SToby Isaac 168720cf1dd8SToby Isaac Input Parameters: 168820cf1dd8SToby Isaac + sp - the PetscDualSpace object 168920cf1dd8SToby Isaac - point - the point (in the dual space's DM) for which the subspace is desired 169020cf1dd8SToby Isaac 169120cf1dd8SToby Isaac Output Parameters: 169220cf1dd8SToby Isaac bdsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the 169320cf1dd8SToby Isaac point, which will be of lesser dimension if height > 0. 169420cf1dd8SToby Isaac 169520cf1dd8SToby Isaac Level: advanced 169620cf1dd8SToby Isaac 169720cf1dd8SToby Isaac .seealso: PetscDualSpace 169820cf1dd8SToby Isaac @*/ 169920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp) 170020cf1dd8SToby Isaac { 1701b4457527SToby Isaac PetscInt pStart = 0, pEnd = 0, cStart, cEnd; 1702b4457527SToby Isaac DM dm; 170320cf1dd8SToby Isaac PetscErrorCode ierr; 170420cf1dd8SToby Isaac 170520cf1dd8SToby Isaac PetscFunctionBegin; 170620cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1707064a246eSJacob Faibussowitsch PetscValidPointer(bdsp,3); 170820cf1dd8SToby Isaac *bdsp = NULL; 1709b4457527SToby Isaac dm = sp->dm; 1710b4457527SToby Isaac ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 17112c71b3e2SJacob Faibussowitsch PetscCheckFalse(point < pStart || point > pEnd,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point"); 1712b4457527SToby Isaac ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 1713b4457527SToby Isaac if (point == cStart && cEnd == cStart + 1) { /* the dual space is only equivalent to the dual space on a cell if the reference mesh has just one cell */ 1714b4457527SToby Isaac *bdsp = sp; 1715b4457527SToby Isaac PetscFunctionReturn(0); 1716b4457527SToby Isaac } 1717b4457527SToby Isaac if (!sp->pointSpaces) { 1718b4457527SToby Isaac PetscInt p; 1719b4457527SToby Isaac ierr = PetscCalloc1(pEnd - pStart, &(sp->pointSpaces));CHKERRQ(ierr); 172020cf1dd8SToby Isaac 1721b4457527SToby Isaac for (p = 0; p < pEnd - pStart; p++) { 1722b4457527SToby Isaac if (p + pStart == cStart && cEnd == cStart + 1) continue; 1723b4457527SToby Isaac if (sp->ops->createpointsubspace) {ierr = (*sp->ops->createpointsubspace)(sp,p+pStart,&(sp->pointSpaces[p]));CHKERRQ(ierr);} 1724b4457527SToby Isaac else if (sp->heightSpaces || sp->ops->createheightsubspace) { 1725b4457527SToby Isaac PetscInt dim, depth, height; 1726b4457527SToby Isaac DMLabel label; 1727b4457527SToby Isaac 172820cf1dd8SToby Isaac ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr); 172920cf1dd8SToby Isaac ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr); 1730b4457527SToby Isaac ierr = DMLabelGetValue(label,p+pStart,&depth);CHKERRQ(ierr); 173120cf1dd8SToby Isaac height = dim - depth; 1732b4457527SToby Isaac ierr = PetscDualSpaceGetHeightSubspace(sp, height, &(sp->pointSpaces[p]));CHKERRQ(ierr); 1733b4457527SToby Isaac ierr = PetscObjectReference((PetscObject)sp->pointSpaces[p]);CHKERRQ(ierr); 173420cf1dd8SToby Isaac } 1735b4457527SToby Isaac } 1736b4457527SToby Isaac } 1737b4457527SToby Isaac *bdsp = sp->pointSpaces[point - pStart]; 173820cf1dd8SToby Isaac PetscFunctionReturn(0); 173920cf1dd8SToby Isaac } 174020cf1dd8SToby Isaac 17416f905325SMatthew G. Knepley /*@C 17426f905325SMatthew G. Knepley PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis 17436f905325SMatthew G. Knepley 17446f905325SMatthew G. Knepley Not collective 17456f905325SMatthew G. Knepley 17466f905325SMatthew G. Knepley Input Parameter: 17476f905325SMatthew G. Knepley . sp - the PetscDualSpace object 17486f905325SMatthew G. Knepley 17496f905325SMatthew G. Knepley Output Parameters: 1750b4457527SToby Isaac + perms - Permutations of the interior degrees of freedom, parameterized by the point orientation 1751b4457527SToby Isaac - flips - Sign reversal of the interior degrees of freedom, parameterized by the point orientation 17526f905325SMatthew G. Knepley 17536f905325SMatthew G. Knepley Note: The permutation and flip arrays are organized in the following way 17546f905325SMatthew G. Knepley $ perms[p][ornt][dof # on point] = new local dof # 17556f905325SMatthew G. Knepley $ flips[p][ornt][dof # on point] = reversal or not 17566f905325SMatthew G. Knepley 17576f905325SMatthew G. Knepley Level: developer 17586f905325SMatthew G. Knepley 17596f905325SMatthew G. Knepley @*/ 17606f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) 17616f905325SMatthew G. Knepley { 17626f905325SMatthew G. Knepley PetscErrorCode ierr; 17636f905325SMatthew G. Knepley 17646f905325SMatthew G. Knepley PetscFunctionBegin; 17656f905325SMatthew G. Knepley PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1); 17666f905325SMatthew G. Knepley if (perms) {PetscValidPointer(perms,2); *perms = NULL;} 17676f905325SMatthew G. Knepley if (flips) {PetscValidPointer(flips,3); *flips = NULL;} 17686f905325SMatthew G. Knepley if (sp->ops->getsymmetries) {ierr = (sp->ops->getsymmetries)(sp,perms,flips);CHKERRQ(ierr);} 17696f905325SMatthew G. Knepley PetscFunctionReturn(0); 17706f905325SMatthew G. Knepley } 17714bee2e38SMatthew G. Knepley 17724bee2e38SMatthew G. Knepley /*@ 1773b4457527SToby Isaac PetscDualSpaceGetFormDegree - Get the form degree k for the k-form the describes the pushforwards/pullbacks of this 1774b4457527SToby Isaac dual space's functionals. 1775b4457527SToby Isaac 1776b4457527SToby Isaac Input Parameter: 1777b4457527SToby Isaac . dsp - The PetscDualSpace 1778b4457527SToby Isaac 1779b4457527SToby Isaac Output Parameter: 1780b4457527SToby Isaac . k - The *signed* degree k of the k. If k >= 0, this means that the degrees of freedom are k-forms, and are stored 1781b4457527SToby Isaac in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms. So for example, 1782b4457527SToby Isaac the 1-form basis in 3-D is (dx, dy, dz), and the 2-form basis in 3-D is (dx wedge dy, dx wedge dz, dy wedge dz). 1783b4457527SToby Isaac If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the 1784b4457527SToby Isaac Hodge star map. So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms 1785b4457527SToby Isaac but are stored as 1-forms. 1786b4457527SToby Isaac 1787b4457527SToby Isaac Level: developer 1788b4457527SToby Isaac 1789b4457527SToby Isaac .seealso: PetscDTAltV, PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType 1790b4457527SToby Isaac @*/ 1791b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace dsp, PetscInt *k) 1792b4457527SToby Isaac { 1793b4457527SToby Isaac PetscFunctionBeginHot; 1794b4457527SToby Isaac PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1795b4457527SToby Isaac PetscValidPointer(k, 2); 1796b4457527SToby Isaac *k = dsp->k; 1797b4457527SToby Isaac PetscFunctionReturn(0); 1798b4457527SToby Isaac } 1799b4457527SToby Isaac 1800b4457527SToby Isaac /*@ 1801b4457527SToby Isaac PetscDualSpaceSetFormDegree - Set the form degree k for the k-form the describes the pushforwards/pullbacks of this 1802b4457527SToby Isaac dual space's functionals. 1803b4457527SToby Isaac 1804d8d19677SJose E. Roman Input Parameters: 1805b4457527SToby Isaac + dsp - The PetscDualSpace 1806b4457527SToby Isaac - k - The *signed* degree k of the k. If k >= 0, this means that the degrees of freedom are k-forms, and are stored 1807b4457527SToby Isaac in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms. So for example, 1808b4457527SToby Isaac the 1-form basis in 3-D is (dx, dy, dz), and the 2-form basis in 3-D is (dx wedge dy, dx wedge dz, dy wedge dz). 1809b4457527SToby Isaac If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the 1810b4457527SToby Isaac Hodge star map. So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms 1811b4457527SToby Isaac but are stored as 1-forms. 1812b4457527SToby Isaac 1813b4457527SToby Isaac Level: developer 1814b4457527SToby Isaac 1815b4457527SToby Isaac .seealso: PetscDTAltV, PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType 1816b4457527SToby Isaac @*/ 1817b4457527SToby Isaac PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace dsp, PetscInt k) 1818b4457527SToby Isaac { 1819b4457527SToby Isaac PetscInt dim; 1820b4457527SToby Isaac 1821b4457527SToby Isaac PetscFunctionBeginHot; 1822b4457527SToby Isaac PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 18232c71b3e2SJacob Faibussowitsch PetscCheckFalse(dsp->setupcalled,PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up"); 1824b4457527SToby Isaac dim = dsp->dm->dim; 18252c71b3e2SJacob Faibussowitsch PetscCheckFalse(k < -dim || k > dim,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported %D-form on %D-dimensional reference cell", PetscAbsInt(k), dim); 1826b4457527SToby Isaac dsp->k = k; 1827b4457527SToby Isaac PetscFunctionReturn(0); 1828b4457527SToby Isaac } 1829b4457527SToby Isaac 1830b4457527SToby Isaac /*@ 18314bee2e38SMatthew G. Knepley PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space 18324bee2e38SMatthew G. Knepley 18334bee2e38SMatthew G. Knepley Input Parameter: 18344bee2e38SMatthew G. Knepley . dsp - The PetscDualSpace 18354bee2e38SMatthew G. Knepley 18364bee2e38SMatthew G. Knepley Output Parameter: 18374bee2e38SMatthew G. Knepley . k - The simplex dimension 18384bee2e38SMatthew G. Knepley 1839a4ce7ad1SMatthew G. Knepley Level: developer 18404bee2e38SMatthew G. Knepley 18414bee2e38SMatthew G. Knepley Note: Currently supported values are 18424bee2e38SMatthew G. Knepley $ 0: These are H_1 methods that only transform coordinates 18434bee2e38SMatthew G. Knepley $ 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM) 18444bee2e38SMatthew G. Knepley $ 2: These are the same as 1 18454bee2e38SMatthew G. Knepley $ 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM) 18464bee2e38SMatthew G. Knepley 18474bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType 18484bee2e38SMatthew G. Knepley @*/ 18494bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k) 18504bee2e38SMatthew G. Knepley { 1851b4457527SToby Isaac PetscInt dim; 1852b4457527SToby Isaac 18534bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 18544bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 18554bee2e38SMatthew G. Knepley PetscValidPointer(k, 2); 1856b4457527SToby Isaac dim = dsp->dm->dim; 1857b4457527SToby Isaac if (!dsp->k) *k = IDENTITY_TRANSFORM; 1858b4457527SToby Isaac else if (dsp->k == 1) *k = COVARIANT_PIOLA_TRANSFORM; 1859b4457527SToby Isaac else if (dsp->k == -(dim - 1)) *k = CONTRAVARIANT_PIOLA_TRANSFORM; 1860b4457527SToby Isaac else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported transformation"); 18614bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 18624bee2e38SMatthew G. Knepley } 18634bee2e38SMatthew G. Knepley 18644bee2e38SMatthew G. Knepley /*@C 18654bee2e38SMatthew G. Knepley PetscDualSpaceTransform - Transform the function values 18664bee2e38SMatthew G. Knepley 18674bee2e38SMatthew G. Knepley Input Parameters: 18684bee2e38SMatthew G. Knepley + dsp - The PetscDualSpace 18694bee2e38SMatthew G. Knepley . trans - The type of transform 18704bee2e38SMatthew G. Knepley . isInverse - Flag to invert the transform 18714bee2e38SMatthew G. Knepley . fegeom - The cell geometry 18724bee2e38SMatthew G. Knepley . Nv - The number of function samples 18734bee2e38SMatthew G. Knepley . Nc - The number of function components 18744bee2e38SMatthew G. Knepley - vals - The function values 18754bee2e38SMatthew G. Knepley 18764bee2e38SMatthew G. Knepley Output Parameter: 18774bee2e38SMatthew G. Knepley . vals - The transformed function values 18784bee2e38SMatthew G. Knepley 1879a4ce7ad1SMatthew G. Knepley Level: intermediate 18804bee2e38SMatthew G. Knepley 1881f9244615SMatthew G. Knepley Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 18822edcad52SToby Isaac 1883f9244615SMatthew G. Knepley .seealso: PetscDualSpaceTransformGradient(), PetscDualSpaceTransformHessian(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType 18844bee2e38SMatthew G. Knepley @*/ 18854bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[]) 18864bee2e38SMatthew G. Knepley { 1887b4457527SToby Isaac PetscReal Jstar[9] = {0}; 1888b4457527SToby Isaac PetscInt dim, v, c, Nk; 1889b4457527SToby Isaac PetscErrorCode ierr; 18904bee2e38SMatthew G. Knepley 18914bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 18924bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 18934bee2e38SMatthew G. Knepley PetscValidPointer(fegeom, 4); 18944bee2e38SMatthew G. Knepley PetscValidPointer(vals, 7); 1895b4457527SToby Isaac /* TODO: not handling dimEmbed != dim right now */ 18962ae266adSMatthew G. Knepley dim = dsp->dm->dim; 1897b4457527SToby Isaac /* No change needed for 0-forms */ 1898b4457527SToby Isaac if (!dsp->k) PetscFunctionReturn(0); 1899b4457527SToby Isaac ierr = PetscDTBinomialInt(dim, PetscAbsInt(dsp->k), &Nk);CHKERRQ(ierr); 1900b4457527SToby Isaac /* TODO: use fegeom->isAffine */ 1901b4457527SToby Isaac ierr = PetscDTAltVPullbackMatrix(dim, dim, isInverse ? fegeom->J : fegeom->invJ, dsp->k, Jstar);CHKERRQ(ierr); 19024bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 1903b4457527SToby Isaac switch (Nk) { 1904b4457527SToby Isaac case 1: 1905b4457527SToby Isaac for (c = 0; c < Nc; c++) vals[v*Nc + c] *= Jstar[0]; 19064bee2e38SMatthew G. Knepley break; 1907b4457527SToby Isaac case 2: 1908b4457527SToby Isaac for (c = 0; c < Nc; c += 2) DMPlex_Mult2DReal_Internal(Jstar, 1, &vals[v*Nc + c], &vals[v*Nc + c]); 19094bee2e38SMatthew G. Knepley break; 1910b4457527SToby Isaac case 3: 1911b4457527SToby Isaac for (c = 0; c < Nc; c += 3) DMPlex_Mult3DReal_Internal(Jstar, 1, &vals[v*Nc + c], &vals[v*Nc + c]); 1912b4457527SToby Isaac break; 191398921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported form size %D for transformation", Nk); 1914b4457527SToby Isaac } 19154bee2e38SMatthew G. Knepley } 19164bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 19174bee2e38SMatthew G. Knepley } 1918b4457527SToby Isaac 19194bee2e38SMatthew G. Knepley /*@C 19204bee2e38SMatthew G. Knepley PetscDualSpaceTransformGradient - Transform the function gradient values 19214bee2e38SMatthew G. Knepley 19224bee2e38SMatthew G. Knepley Input Parameters: 19234bee2e38SMatthew G. Knepley + dsp - The PetscDualSpace 19244bee2e38SMatthew G. Knepley . trans - The type of transform 19254bee2e38SMatthew G. Knepley . isInverse - Flag to invert the transform 19264bee2e38SMatthew G. Knepley . fegeom - The cell geometry 19274bee2e38SMatthew G. Knepley . Nv - The number of function gradient samples 19284bee2e38SMatthew G. Knepley . Nc - The number of function components 19294bee2e38SMatthew G. Knepley - vals - The function gradient values 19304bee2e38SMatthew G. Knepley 19314bee2e38SMatthew G. Knepley Output Parameter: 1932f9244615SMatthew G. Knepley . vals - The transformed function gradient values 19334bee2e38SMatthew G. Knepley 1934a4ce7ad1SMatthew G. Knepley Level: intermediate 19354bee2e38SMatthew G. Knepley 1936f9244615SMatthew G. Knepley Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 19372edcad52SToby Isaac 1938625e0966SMatthew G. Knepley .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType 19394bee2e38SMatthew G. Knepley @*/ 19404bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[]) 19414bee2e38SMatthew G. Knepley { 194227f02ce8SMatthew G. Knepley const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed; 194327f02ce8SMatthew G. Knepley PetscInt v, c, d; 19444bee2e38SMatthew G. Knepley 19454bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 19464bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 19474bee2e38SMatthew G. Knepley PetscValidPointer(fegeom, 4); 19484bee2e38SMatthew G. Knepley PetscValidPointer(vals, 7); 194927f02ce8SMatthew G. Knepley #ifdef PETSC_USE_DEBUG 19502c71b3e2SJacob Faibussowitsch PetscCheckFalse(dE <= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %D", dE); 195127f02ce8SMatthew G. Knepley #endif 19524bee2e38SMatthew G. Knepley /* Transform gradient */ 195327f02ce8SMatthew G. Knepley if (dim == dE) { 19544bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 19554bee2e38SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 19564bee2e38SMatthew G. Knepley switch (dim) 19574bee2e38SMatthew G. Knepley { 1958100a78e1SStefano Zampini case 1: vals[(v*Nc+c)*dim] *= fegeom->invJ[0];break; 19596142fa51SMatthew G. Knepley case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break; 19606142fa51SMatthew G. Knepley case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break; 196198921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 19624bee2e38SMatthew G. Knepley } 19634bee2e38SMatthew G. Knepley } 19644bee2e38SMatthew G. Knepley } 196527f02ce8SMatthew G. Knepley } else { 196627f02ce8SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 196727f02ce8SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 196827f02ce8SMatthew G. Knepley DMPlex_MultTransposeReal_Internal(fegeom->invJ, dim, dE, 1, &vals[(v*Nc+c)*dE], &vals[(v*Nc+c)*dE]); 196927f02ce8SMatthew G. Knepley } 197027f02ce8SMatthew G. Knepley } 197127f02ce8SMatthew G. Knepley } 19724bee2e38SMatthew G. Knepley /* Assume its a vector, otherwise assume its a bunch of scalars */ 19734bee2e38SMatthew G. Knepley if (Nc == 1 || Nc != dim) PetscFunctionReturn(0); 19744bee2e38SMatthew G. Knepley switch (trans) { 19754bee2e38SMatthew G. Knepley case IDENTITY_TRANSFORM: break; 19764bee2e38SMatthew G. Knepley case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */ 19774bee2e38SMatthew G. Knepley if (isInverse) { 19784bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 19794bee2e38SMatthew G. Knepley for (d = 0; d < dim; ++d) { 19804bee2e38SMatthew G. Knepley switch (dim) 19814bee2e38SMatthew G. Knepley { 19826142fa51SMatthew G. Knepley case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 19836142fa51SMatthew G. Knepley case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 198498921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 19854bee2e38SMatthew G. Knepley } 19864bee2e38SMatthew G. Knepley } 19874bee2e38SMatthew G. Knepley } 19884bee2e38SMatthew G. Knepley } else { 19894bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 19904bee2e38SMatthew G. Knepley for (d = 0; d < dim; ++d) { 19914bee2e38SMatthew G. Knepley switch (dim) 19924bee2e38SMatthew G. Knepley { 19936142fa51SMatthew G. Knepley case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 19946142fa51SMatthew G. Knepley case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 199598921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 19964bee2e38SMatthew G. Knepley } 19974bee2e38SMatthew G. Knepley } 19984bee2e38SMatthew G. Knepley } 19994bee2e38SMatthew G. Knepley } 20004bee2e38SMatthew G. Knepley break; 20014bee2e38SMatthew G. Knepley case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */ 20024bee2e38SMatthew G. Knepley if (isInverse) { 20034bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 20044bee2e38SMatthew G. Knepley for (d = 0; d < dim; ++d) { 20054bee2e38SMatthew G. Knepley switch (dim) 20064bee2e38SMatthew G. Knepley { 20076142fa51SMatthew G. Knepley case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 20086142fa51SMatthew G. Knepley case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 200998921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 20104bee2e38SMatthew G. Knepley } 20114bee2e38SMatthew G. Knepley for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] *= fegeom->detJ[0]; 20124bee2e38SMatthew G. Knepley } 20134bee2e38SMatthew G. Knepley } 20144bee2e38SMatthew G. Knepley } else { 20154bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 20164bee2e38SMatthew G. Knepley for (d = 0; d < dim; ++d) { 20174bee2e38SMatthew G. Knepley switch (dim) 20184bee2e38SMatthew G. Knepley { 20196142fa51SMatthew G. Knepley case 2: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 20206142fa51SMatthew G. Knepley case 3: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 202198921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 20224bee2e38SMatthew G. Knepley } 20234bee2e38SMatthew G. Knepley for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] /= fegeom->detJ[0]; 20244bee2e38SMatthew G. Knepley } 20254bee2e38SMatthew G. Knepley } 20264bee2e38SMatthew G. Knepley } 20274bee2e38SMatthew G. Knepley break; 20284bee2e38SMatthew G. Knepley } 20294bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 20304bee2e38SMatthew G. Knepley } 20314bee2e38SMatthew G. Knepley 20324bee2e38SMatthew G. Knepley /*@C 2033f9244615SMatthew G. Knepley PetscDualSpaceTransformHessian - Transform the function Hessian values 2034f9244615SMatthew G. Knepley 2035f9244615SMatthew G. Knepley Input Parameters: 2036f9244615SMatthew G. Knepley + dsp - The PetscDualSpace 2037f9244615SMatthew G. Knepley . trans - The type of transform 2038f9244615SMatthew G. Knepley . isInverse - Flag to invert the transform 2039f9244615SMatthew G. Knepley . fegeom - The cell geometry 2040f9244615SMatthew G. Knepley . Nv - The number of function Hessian samples 2041f9244615SMatthew G. Knepley . Nc - The number of function components 2042f9244615SMatthew G. Knepley - vals - The function gradient values 2043f9244615SMatthew G. Knepley 2044f9244615SMatthew G. Knepley Output Parameter: 2045f9244615SMatthew G. Knepley . vals - The transformed function Hessian values 2046f9244615SMatthew G. Knepley 2047f9244615SMatthew G. Knepley Level: intermediate 2048f9244615SMatthew G. Knepley 2049f9244615SMatthew G. Knepley Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 2050f9244615SMatthew G. Knepley 2051f9244615SMatthew G. Knepley .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType 2052f9244615SMatthew G. Knepley @*/ 2053f9244615SMatthew G. Knepley PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[]) 2054f9244615SMatthew G. Knepley { 2055f9244615SMatthew G. Knepley const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed; 2056f9244615SMatthew G. Knepley PetscInt v, c; 2057f9244615SMatthew G. Knepley 2058f9244615SMatthew G. Knepley PetscFunctionBeginHot; 2059f9244615SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 2060f9244615SMatthew G. Knepley PetscValidPointer(fegeom, 4); 2061f9244615SMatthew G. Knepley PetscValidPointer(vals, 7); 2062f9244615SMatthew G. Knepley #ifdef PETSC_USE_DEBUG 20632c71b3e2SJacob Faibussowitsch PetscCheckFalse(dE <= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %D", dE); 2064f9244615SMatthew G. Knepley #endif 2065f9244615SMatthew G. Knepley /* Transform Hessian: J^{-T}_{ik} J^{-T}_{jl} H(f)_{kl} = J^{-T}_{ik} H(f)_{kl} J^{-1}_{lj} */ 2066f9244615SMatthew G. Knepley if (dim == dE) { 2067f9244615SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 2068f9244615SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 2069f9244615SMatthew G. Knepley switch (dim) 2070f9244615SMatthew G. Knepley { 2071f9244615SMatthew G. Knepley case 1: vals[(v*Nc+c)*dim*dim] *= PetscSqr(fegeom->invJ[0]);break; 2072f9244615SMatthew G. Knepley case 2: DMPlex_PTAP2DReal_Internal(fegeom->invJ, &vals[(v*Nc+c)*dim*dim], &vals[(v*Nc+c)*dim*dim]);break; 2073f9244615SMatthew G. Knepley case 3: DMPlex_PTAP3DReal_Internal(fegeom->invJ, &vals[(v*Nc+c)*dim*dim], &vals[(v*Nc+c)*dim*dim]);break; 207498921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 2075f9244615SMatthew G. Knepley } 2076f9244615SMatthew G. Knepley } 2077f9244615SMatthew G. Knepley } 2078f9244615SMatthew G. Knepley } else { 2079f9244615SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 2080f9244615SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 2081f9244615SMatthew G. Knepley DMPlex_PTAPReal_Internal(fegeom->invJ, dim, dE, &vals[(v*Nc+c)*dE*dE], &vals[(v*Nc+c)*dE*dE]); 2082f9244615SMatthew G. Knepley } 2083f9244615SMatthew G. Knepley } 2084f9244615SMatthew G. Knepley } 2085f9244615SMatthew G. Knepley /* Assume its a vector, otherwise assume its a bunch of scalars */ 2086f9244615SMatthew G. Knepley if (Nc == 1 || Nc != dim) PetscFunctionReturn(0); 2087f9244615SMatthew G. Knepley switch (trans) { 2088f9244615SMatthew G. Knepley case IDENTITY_TRANSFORM: break; 2089f9244615SMatthew G. Knepley case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */ 2090f9244615SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported"); 2091f9244615SMatthew G. Knepley case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */ 2092f9244615SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported"); 2093f9244615SMatthew G. Knepley } 2094f9244615SMatthew G. Knepley PetscFunctionReturn(0); 2095f9244615SMatthew G. Knepley } 2096f9244615SMatthew G. Knepley 2097f9244615SMatthew G. Knepley /*@C 20984bee2e38SMatthew G. Knepley PetscDualSpacePullback - Transform the given functional so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method. 20994bee2e38SMatthew G. Knepley 21004bee2e38SMatthew G. Knepley Input Parameters: 21014bee2e38SMatthew G. Knepley + dsp - The PetscDualSpace 21024bee2e38SMatthew G. Knepley . fegeom - The geometry for this cell 21034bee2e38SMatthew G. Knepley . Nq - The number of function samples 21044bee2e38SMatthew G. Knepley . Nc - The number of function components 21054bee2e38SMatthew G. Knepley - pointEval - The function values 21064bee2e38SMatthew G. Knepley 21074bee2e38SMatthew G. Knepley Output Parameter: 21084bee2e38SMatthew G. Knepley . pointEval - The transformed function values 21094bee2e38SMatthew G. Knepley 21104bee2e38SMatthew G. Knepley Level: advanced 21114bee2e38SMatthew G. Knepley 21124bee2e38SMatthew G. Knepley Note: Functions transform in a complementary way (pushforward) to functionals, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex. 21134bee2e38SMatthew G. Knepley 21142edcad52SToby Isaac Note: This only handles tranformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 21152edcad52SToby Isaac 21164bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 21174bee2e38SMatthew G. Knepley @*/ 21182edcad52SToby Isaac PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 21194bee2e38SMatthew G. Knepley { 21204bee2e38SMatthew G. Knepley PetscDualSpaceTransformType trans; 2121b4457527SToby Isaac PetscInt k; 21224bee2e38SMatthew G. Knepley PetscErrorCode ierr; 21234bee2e38SMatthew G. Knepley 21244bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 21254bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 21264bee2e38SMatthew G. Knepley PetscValidPointer(fegeom, 2); 21272edcad52SToby Isaac PetscValidPointer(pointEval, 5); 21284bee2e38SMatthew G. Knepley /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 21294bee2e38SMatthew G. Knepley This determines their transformation properties. */ 2130b4457527SToby Isaac ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr); 2131b4457527SToby Isaac switch (k) 21324bee2e38SMatthew G. Knepley { 21334bee2e38SMatthew G. Knepley case 0: /* H^1 point evaluations */ 21344bee2e38SMatthew G. Knepley trans = IDENTITY_TRANSFORM;break; 21354bee2e38SMatthew G. Knepley case 1: /* Hcurl preserves tangential edge traces */ 21364bee2e38SMatthew G. Knepley trans = COVARIANT_PIOLA_TRANSFORM;break; 2137b4457527SToby Isaac case 2: 21384bee2e38SMatthew G. Knepley case 3: /* Hdiv preserve normal traces */ 21394bee2e38SMatthew G. Knepley trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 214098921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k); 21414bee2e38SMatthew G. Knepley } 21422edcad52SToby Isaac ierr = PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr); 21434bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 21444bee2e38SMatthew G. Knepley } 21454bee2e38SMatthew G. Knepley 21464bee2e38SMatthew G. Knepley /*@C 21474bee2e38SMatthew G. Knepley PetscDualSpacePushforward - Transform the given function so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method. 21484bee2e38SMatthew G. Knepley 21494bee2e38SMatthew G. Knepley Input Parameters: 21504bee2e38SMatthew G. Knepley + dsp - The PetscDualSpace 21514bee2e38SMatthew G. Knepley . fegeom - The geometry for this cell 21524bee2e38SMatthew G. Knepley . Nq - The number of function samples 21534bee2e38SMatthew G. Knepley . Nc - The number of function components 21544bee2e38SMatthew G. Knepley - pointEval - The function values 21554bee2e38SMatthew G. Knepley 21564bee2e38SMatthew G. Knepley Output Parameter: 21574bee2e38SMatthew G. Knepley . pointEval - The transformed function values 21584bee2e38SMatthew G. Knepley 21594bee2e38SMatthew G. Knepley Level: advanced 21604bee2e38SMatthew G. Knepley 21614bee2e38SMatthew G. Knepley Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex. 21624bee2e38SMatthew G. Knepley 2163f9244615SMatthew G. Knepley Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 21642edcad52SToby Isaac 21654bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 21664bee2e38SMatthew G. Knepley @*/ 21672edcad52SToby Isaac PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 21684bee2e38SMatthew G. Knepley { 21694bee2e38SMatthew G. Knepley PetscDualSpaceTransformType trans; 2170b4457527SToby Isaac PetscInt k; 21714bee2e38SMatthew G. Knepley PetscErrorCode ierr; 21724bee2e38SMatthew G. Knepley 21734bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 21744bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 21754bee2e38SMatthew G. Knepley PetscValidPointer(fegeom, 2); 21762edcad52SToby Isaac PetscValidPointer(pointEval, 5); 21774bee2e38SMatthew G. Knepley /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 21784bee2e38SMatthew G. Knepley This determines their transformation properties. */ 2179b4457527SToby Isaac ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr); 2180b4457527SToby Isaac switch (k) 21814bee2e38SMatthew G. Knepley { 21824bee2e38SMatthew G. Knepley case 0: /* H^1 point evaluations */ 21834bee2e38SMatthew G. Knepley trans = IDENTITY_TRANSFORM;break; 21844bee2e38SMatthew G. Knepley case 1: /* Hcurl preserves tangential edge traces */ 21854bee2e38SMatthew G. Knepley trans = COVARIANT_PIOLA_TRANSFORM;break; 2186b4457527SToby Isaac case 2: 21874bee2e38SMatthew G. Knepley case 3: /* Hdiv preserve normal traces */ 21884bee2e38SMatthew G. Knepley trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 218998921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k); 21904bee2e38SMatthew G. Knepley } 21912edcad52SToby Isaac ierr = PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr); 21924bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 21934bee2e38SMatthew G. Knepley } 21944bee2e38SMatthew G. Knepley 21954bee2e38SMatthew G. Knepley /*@C 21964bee2e38SMatthew G. Knepley PetscDualSpacePushforwardGradient - Transform the given function gradient so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method. 21974bee2e38SMatthew G. Knepley 21984bee2e38SMatthew G. Knepley Input Parameters: 21994bee2e38SMatthew G. Knepley + dsp - The PetscDualSpace 22004bee2e38SMatthew G. Knepley . fegeom - The geometry for this cell 22014bee2e38SMatthew G. Knepley . Nq - The number of function gradient samples 22024bee2e38SMatthew G. Knepley . Nc - The number of function components 22034bee2e38SMatthew G. Knepley - pointEval - The function gradient values 22044bee2e38SMatthew G. Knepley 22054bee2e38SMatthew G. Knepley Output Parameter: 22064bee2e38SMatthew G. Knepley . pointEval - The transformed function gradient values 22074bee2e38SMatthew G. Knepley 22084bee2e38SMatthew G. Knepley Level: advanced 22094bee2e38SMatthew G. Knepley 22104bee2e38SMatthew G. Knepley Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex. 22114bee2e38SMatthew G. Knepley 2212f9244615SMatthew G. Knepley Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 22132edcad52SToby Isaac 22144bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 2215dc0529c6SBarry Smith @*/ 22162edcad52SToby Isaac PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 22174bee2e38SMatthew G. Knepley { 22184bee2e38SMatthew G. Knepley PetscDualSpaceTransformType trans; 2219b4457527SToby Isaac PetscInt k; 22204bee2e38SMatthew G. Knepley PetscErrorCode ierr; 22214bee2e38SMatthew G. Knepley 22224bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 22234bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 22244bee2e38SMatthew G. Knepley PetscValidPointer(fegeom, 2); 22252edcad52SToby Isaac PetscValidPointer(pointEval, 5); 22264bee2e38SMatthew G. Knepley /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 22274bee2e38SMatthew G. Knepley This determines their transformation properties. */ 2228b4457527SToby Isaac ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr); 2229b4457527SToby Isaac switch (k) 22304bee2e38SMatthew G. Knepley { 22314bee2e38SMatthew G. Knepley case 0: /* H^1 point evaluations */ 22324bee2e38SMatthew G. Knepley trans = IDENTITY_TRANSFORM;break; 22334bee2e38SMatthew G. Knepley case 1: /* Hcurl preserves tangential edge traces */ 22344bee2e38SMatthew G. Knepley trans = COVARIANT_PIOLA_TRANSFORM;break; 2235b4457527SToby Isaac case 2: 22364bee2e38SMatthew G. Knepley case 3: /* Hdiv preserve normal traces */ 22374bee2e38SMatthew G. Knepley trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 223898921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k); 22394bee2e38SMatthew G. Knepley } 22402edcad52SToby Isaac ierr = PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr); 22414bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 22424bee2e38SMatthew G. Knepley } 2243f9244615SMatthew G. Knepley 2244f9244615SMatthew G. Knepley /*@C 2245f9244615SMatthew G. Knepley PetscDualSpacePushforwardHessian - Transform the given function Hessian so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method. 2246f9244615SMatthew G. Knepley 2247f9244615SMatthew G. Knepley Input Parameters: 2248f9244615SMatthew G. Knepley + dsp - The PetscDualSpace 2249f9244615SMatthew G. Knepley . fegeom - The geometry for this cell 2250f9244615SMatthew G. Knepley . Nq - The number of function Hessian samples 2251f9244615SMatthew G. Knepley . Nc - The number of function components 2252f9244615SMatthew G. Knepley - pointEval - The function gradient values 2253f9244615SMatthew G. Knepley 2254f9244615SMatthew G. Knepley Output Parameter: 2255f9244615SMatthew G. Knepley . pointEval - The transformed function Hessian values 2256f9244615SMatthew G. Knepley 2257f9244615SMatthew G. Knepley Level: advanced 2258f9244615SMatthew G. Knepley 2259f9244615SMatthew G. Knepley Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex. 2260f9244615SMatthew G. Knepley 2261f9244615SMatthew G. Knepley Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 2262f9244615SMatthew G. Knepley 2263f9244615SMatthew G. Knepley .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 2264f9244615SMatthew G. Knepley @*/ 2265f9244615SMatthew G. Knepley PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 2266f9244615SMatthew G. Knepley { 2267f9244615SMatthew G. Knepley PetscDualSpaceTransformType trans; 2268f9244615SMatthew G. Knepley PetscInt k; 2269f9244615SMatthew G. Knepley PetscErrorCode ierr; 2270f9244615SMatthew G. Knepley 2271f9244615SMatthew G. Knepley PetscFunctionBeginHot; 2272f9244615SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 2273f9244615SMatthew G. Knepley PetscValidPointer(fegeom, 2); 2274f9244615SMatthew G. Knepley PetscValidPointer(pointEval, 5); 2275f9244615SMatthew G. Knepley /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 2276f9244615SMatthew G. Knepley This determines their transformation properties. */ 2277f9244615SMatthew G. Knepley ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr); 2278f9244615SMatthew G. Knepley switch (k) 2279f9244615SMatthew G. Knepley { 2280f9244615SMatthew G. Knepley case 0: /* H^1 point evaluations */ 2281f9244615SMatthew G. Knepley trans = IDENTITY_TRANSFORM;break; 2282f9244615SMatthew G. Knepley case 1: /* Hcurl preserves tangential edge traces */ 2283f9244615SMatthew G. Knepley trans = COVARIANT_PIOLA_TRANSFORM;break; 2284f9244615SMatthew G. Knepley case 2: 2285f9244615SMatthew G. Knepley case 3: /* Hdiv preserve normal traces */ 2286f9244615SMatthew G. Knepley trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 228798921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k); 2288f9244615SMatthew G. Knepley } 2289f9244615SMatthew G. Knepley ierr = PetscDualSpaceTransformHessian(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr); 2290f9244615SMatthew G. Knepley PetscFunctionReturn(0); 2291f9244615SMatthew G. Knepley } 2292