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 620cf1dd8SToby Isaac PetscFunctionList PetscDualSpaceList = NULL; 720cf1dd8SToby Isaac PetscBool PetscDualSpaceRegisterAllCalled = PETSC_FALSE; 820cf1dd8SToby Isaac 955cc6565SMatthew G. Knepley const char *const PetscDualSpaceReferenceCells[] = {"SIMPLEX", "TENSOR", "PetscDualSpaceReferenceCell", "PETSCDUALSPACE_REFCELL_",0}; 1055cc6565SMatthew G. Knepley 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. 146f905325SMatthew 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,0}. 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); 14520cf1dd8SToby Isaac if (!r) SETERRQ1(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); 195221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim);CHKERRQ(ierr); 196221d6281SMatthew G. Knepley if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);} 197221d6281SMatthew G. Knepley ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr); 198221d6281SMatthew G. Knepley if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 199221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 200221d6281SMatthew G. Knepley for (f = 0; f < pdim; ++f) { 201221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f);CHKERRQ(ierr); 202221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 203221d6281SMatthew G. Knepley ierr = PetscQuadratureView(sp->functional[f], v);CHKERRQ(ierr); 204221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 205221d6281SMatthew G. Knepley } 206221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 207221d6281SMatthew G. Knepley } 208221d6281SMatthew G. Knepley ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 209221d6281SMatthew G. Knepley PetscFunctionReturn(0); 210221d6281SMatthew G. Knepley } 211221d6281SMatthew G. Knepley 212fe2efc57SMark /*@C 213fe2efc57SMark PetscDualSpaceViewFromOptions - View from Options 214fe2efc57SMark 215fe2efc57SMark Collective on PetscDualSpace 216fe2efc57SMark 217fe2efc57SMark Input Parameters: 218fe2efc57SMark + A - the PetscDualSpace object 219736c3998SJose E. Roman . obj - Optional object, proivides prefix 220736c3998SJose E. Roman - name - command line option 221fe2efc57SMark 222fe2efc57SMark Level: intermediate 223fe2efc57SMark .seealso: PetscDualSpace, PetscDualSpaceView(), PetscObjectViewFromOptions(), PetscDualSpaceCreate() 224fe2efc57SMark @*/ 225fe2efc57SMark PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace A,PetscObject obj,const char name[]) 226fe2efc57SMark { 227fe2efc57SMark PetscErrorCode ierr; 228fe2efc57SMark 229fe2efc57SMark PetscFunctionBegin; 230fe2efc57SMark PetscValidHeaderSpecific(A,PETSCDUALSPACE_CLASSID,1); 231fe2efc57SMark ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr); 232fe2efc57SMark PetscFunctionReturn(0); 233fe2efc57SMark } 234fe2efc57SMark 23520cf1dd8SToby Isaac /*@ 23620cf1dd8SToby Isaac PetscDualSpaceView - Views a PetscDualSpace 23720cf1dd8SToby Isaac 238d083f849SBarry Smith Collective on sp 23920cf1dd8SToby Isaac 24020cf1dd8SToby Isaac Input Parameter: 24120cf1dd8SToby Isaac + sp - the PetscDualSpace object to view 24220cf1dd8SToby Isaac - v - the viewer 24320cf1dd8SToby Isaac 244a4ce7ad1SMatthew G. Knepley Level: beginner 24520cf1dd8SToby Isaac 246fe2efc57SMark .seealso PetscDualSpaceDestroy(), PetscDualSpace 24720cf1dd8SToby Isaac @*/ 24820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v) 24920cf1dd8SToby Isaac { 250d9bac1caSLisandro Dalcin PetscBool iascii; 25120cf1dd8SToby Isaac PetscErrorCode ierr; 25220cf1dd8SToby Isaac 25320cf1dd8SToby Isaac PetscFunctionBegin; 25420cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 255d9bac1caSLisandro Dalcin if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2); 25620cf1dd8SToby Isaac if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);} 257d9bac1caSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 258221d6281SMatthew G. Knepley if (iascii) {ierr = PetscDualSpaceView_ASCII(sp, v);CHKERRQ(ierr);} 25920cf1dd8SToby Isaac PetscFunctionReturn(0); 26020cf1dd8SToby Isaac } 26120cf1dd8SToby Isaac 26220cf1dd8SToby Isaac /*@ 26320cf1dd8SToby Isaac PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database 26420cf1dd8SToby Isaac 265d083f849SBarry Smith Collective on sp 26620cf1dd8SToby Isaac 26720cf1dd8SToby Isaac Input Parameter: 26820cf1dd8SToby Isaac . sp - the PetscDualSpace object to set options for 26920cf1dd8SToby Isaac 27020cf1dd8SToby Isaac Options Database: 2717be5e748SToby Isaac . -petscspace_degree the approximation order of the space 27220cf1dd8SToby Isaac 273a4ce7ad1SMatthew G. Knepley Level: intermediate 27420cf1dd8SToby Isaac 275fe2efc57SMark .seealso PetscDualSpaceView(), PetscDualSpace, PetscObjectSetFromOptions() 27620cf1dd8SToby Isaac @*/ 27720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp) 27820cf1dd8SToby Isaac { 279063ee4adSMatthew G. Knepley PetscDualSpaceReferenceCell refCell = PETSCDUALSPACE_REFCELL_SIMPLEX; 280063ee4adSMatthew G. Knepley PetscInt refDim = 0; 281063ee4adSMatthew G. Knepley PetscBool flg; 28220cf1dd8SToby Isaac const char *defaultType; 28320cf1dd8SToby Isaac char name[256]; 28420cf1dd8SToby Isaac PetscErrorCode ierr; 28520cf1dd8SToby Isaac 28620cf1dd8SToby Isaac PetscFunctionBegin; 28720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 28820cf1dd8SToby Isaac if (!((PetscObject) sp)->type_name) { 28920cf1dd8SToby Isaac defaultType = PETSCDUALSPACELAGRANGE; 29020cf1dd8SToby Isaac } else { 29120cf1dd8SToby Isaac defaultType = ((PetscObject) sp)->type_name; 29220cf1dd8SToby Isaac } 29320cf1dd8SToby Isaac if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 29420cf1dd8SToby Isaac 29520cf1dd8SToby Isaac ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 29620cf1dd8SToby Isaac ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr); 29720cf1dd8SToby Isaac if (flg) { 29820cf1dd8SToby Isaac ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr); 29920cf1dd8SToby Isaac } else if (!((PetscObject) sp)->type_name) { 30020cf1dd8SToby Isaac ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr); 30120cf1dd8SToby Isaac } 3025a856986SBarry Smith ierr = PetscOptionsBoundedInt("-petscdualspace_degree", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL,0);CHKERRQ(ierr); 3035a856986SBarry Smith ierr = PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL,1);CHKERRQ(ierr); 30420cf1dd8SToby Isaac if (sp->ops->setfromoptions) { 30520cf1dd8SToby Isaac ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr); 30620cf1dd8SToby Isaac } 3075a856986SBarry Smith ierr = PetscOptionsBoundedInt("-petscdualspace_refdim", "The spatial dimension of the reference cell", "PetscDualSpaceSetReferenceCell", refDim, &refDim, NULL,0);CHKERRQ(ierr); 308063ee4adSMatthew G. Knepley ierr = PetscOptionsEnum("-petscdualspace_refcell", "Reference cell", "PetscDualSpaceSetReferenceCell", PetscDualSpaceReferenceCells, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);CHKERRQ(ierr); 309063ee4adSMatthew G. Knepley if (flg) { 310063ee4adSMatthew G. Knepley DM K; 311063ee4adSMatthew G. Knepley 312063ee4adSMatthew G. Knepley if (!refDim) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_INCOMP, "Reference cell specified without a dimension. Use -petscdualspace_refdim."); 313063ee4adSMatthew G. Knepley ierr = PetscDualSpaceCreateReferenceCell(sp, refDim, refCell == PETSCDUALSPACE_REFCELL_SIMPLEX ? PETSC_TRUE : PETSC_FALSE, &K);CHKERRQ(ierr); 314063ee4adSMatthew G. Knepley ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr); 315063ee4adSMatthew G. Knepley ierr = DMDestroy(&K);CHKERRQ(ierr); 316063ee4adSMatthew G. Knepley } 317063ee4adSMatthew G. Knepley 31820cf1dd8SToby Isaac /* process any options handlers added with PetscObjectAddOptionsHandler() */ 31920cf1dd8SToby Isaac ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr); 32020cf1dd8SToby Isaac ierr = PetscOptionsEnd();CHKERRQ(ierr); 321063ee4adSMatthew G. Knepley sp->setfromoptionscalled = PETSC_TRUE; 32220cf1dd8SToby Isaac PetscFunctionReturn(0); 32320cf1dd8SToby Isaac } 32420cf1dd8SToby Isaac 32520cf1dd8SToby Isaac /*@ 32620cf1dd8SToby Isaac PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace 32720cf1dd8SToby Isaac 328d083f849SBarry Smith Collective on sp 32920cf1dd8SToby Isaac 33020cf1dd8SToby Isaac Input Parameter: 33120cf1dd8SToby Isaac . sp - the PetscDualSpace object to setup 33220cf1dd8SToby Isaac 333a4ce7ad1SMatthew G. Knepley Level: intermediate 33420cf1dd8SToby Isaac 335fe2efc57SMark .seealso PetscDualSpaceView(), PetscDualSpaceDestroy(), PetscDualSpace 33620cf1dd8SToby Isaac @*/ 33720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp) 33820cf1dd8SToby Isaac { 33920cf1dd8SToby Isaac PetscErrorCode ierr; 34020cf1dd8SToby Isaac 34120cf1dd8SToby Isaac PetscFunctionBegin; 34220cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 34320cf1dd8SToby Isaac if (sp->setupcalled) PetscFunctionReturn(0); 34420cf1dd8SToby Isaac sp->setupcalled = PETSC_TRUE; 34520cf1dd8SToby Isaac if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);} 346063ee4adSMatthew G. Knepley if (sp->setfromoptionscalled) {ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);} 34720cf1dd8SToby Isaac PetscFunctionReturn(0); 34820cf1dd8SToby Isaac } 34920cf1dd8SToby Isaac 35020cf1dd8SToby Isaac /*@ 35120cf1dd8SToby Isaac PetscDualSpaceDestroy - Destroys a PetscDualSpace object 35220cf1dd8SToby Isaac 353d083f849SBarry Smith Collective on sp 35420cf1dd8SToby Isaac 35520cf1dd8SToby Isaac Input Parameter: 35620cf1dd8SToby Isaac . sp - the PetscDualSpace object to destroy 35720cf1dd8SToby Isaac 358a4ce7ad1SMatthew G. Knepley Level: beginner 35920cf1dd8SToby Isaac 360fe2efc57SMark .seealso PetscDualSpaceView(), PetscDualSpace(), PetscDualSpaceCreate() 36120cf1dd8SToby Isaac @*/ 36220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp) 36320cf1dd8SToby Isaac { 36420cf1dd8SToby Isaac PetscInt dim, f; 36520cf1dd8SToby Isaac PetscErrorCode ierr; 36620cf1dd8SToby Isaac 36720cf1dd8SToby Isaac PetscFunctionBegin; 36820cf1dd8SToby Isaac if (!*sp) PetscFunctionReturn(0); 36920cf1dd8SToby Isaac PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1); 37020cf1dd8SToby Isaac 37120cf1dd8SToby Isaac if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);} 37220cf1dd8SToby Isaac ((PetscObject) (*sp))->refct = 0; 37320cf1dd8SToby Isaac 37420cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr); 37520cf1dd8SToby Isaac for (f = 0; f < dim; ++f) { 37620cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr); 37720cf1dd8SToby Isaac } 37820cf1dd8SToby Isaac ierr = PetscFree((*sp)->functional);CHKERRQ(ierr); 37920cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&(*sp)->allPoints);CHKERRQ(ierr); 38020cf1dd8SToby Isaac ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr); 38120cf1dd8SToby Isaac 38220cf1dd8SToby Isaac if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);} 38320cf1dd8SToby Isaac ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr); 38420cf1dd8SToby Isaac PetscFunctionReturn(0); 38520cf1dd8SToby Isaac } 38620cf1dd8SToby Isaac 38720cf1dd8SToby Isaac /*@ 38820cf1dd8SToby Isaac PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType(). 38920cf1dd8SToby Isaac 390d083f849SBarry Smith Collective 39120cf1dd8SToby Isaac 39220cf1dd8SToby Isaac Input Parameter: 39320cf1dd8SToby Isaac . comm - The communicator for the PetscDualSpace object 39420cf1dd8SToby Isaac 39520cf1dd8SToby Isaac Output Parameter: 39620cf1dd8SToby Isaac . sp - The PetscDualSpace object 39720cf1dd8SToby Isaac 39820cf1dd8SToby Isaac Level: beginner 39920cf1dd8SToby Isaac 40020cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE 40120cf1dd8SToby Isaac @*/ 40220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp) 40320cf1dd8SToby Isaac { 40420cf1dd8SToby Isaac PetscDualSpace s; 40520cf1dd8SToby Isaac PetscErrorCode ierr; 40620cf1dd8SToby Isaac 40720cf1dd8SToby Isaac PetscFunctionBegin; 40820cf1dd8SToby Isaac PetscValidPointer(sp, 2); 40920cf1dd8SToby Isaac ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr); 41020cf1dd8SToby Isaac *sp = NULL; 41120cf1dd8SToby Isaac ierr = PetscFEInitializePackage();CHKERRQ(ierr); 41220cf1dd8SToby Isaac 41320cf1dd8SToby Isaac ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr); 41420cf1dd8SToby Isaac 41520cf1dd8SToby Isaac s->order = 0; 41620cf1dd8SToby Isaac s->Nc = 1; 4174bee2e38SMatthew G. Knepley s->k = 0; 41820cf1dd8SToby Isaac s->setupcalled = PETSC_FALSE; 41920cf1dd8SToby Isaac 42020cf1dd8SToby Isaac *sp = s; 42120cf1dd8SToby Isaac PetscFunctionReturn(0); 42220cf1dd8SToby Isaac } 42320cf1dd8SToby Isaac 42420cf1dd8SToby Isaac /*@ 42520cf1dd8SToby Isaac PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup. 42620cf1dd8SToby Isaac 427d083f849SBarry Smith Collective on sp 42820cf1dd8SToby Isaac 42920cf1dd8SToby Isaac Input Parameter: 43020cf1dd8SToby Isaac . sp - The original PetscDualSpace 43120cf1dd8SToby Isaac 43220cf1dd8SToby Isaac Output Parameter: 43320cf1dd8SToby Isaac . spNew - The duplicate PetscDualSpace 43420cf1dd8SToby Isaac 43520cf1dd8SToby Isaac Level: beginner 43620cf1dd8SToby Isaac 43720cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType() 43820cf1dd8SToby Isaac @*/ 43920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew) 44020cf1dd8SToby Isaac { 44120cf1dd8SToby Isaac PetscErrorCode ierr; 44220cf1dd8SToby Isaac 44320cf1dd8SToby Isaac PetscFunctionBegin; 44420cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 44520cf1dd8SToby Isaac PetscValidPointer(spNew, 2); 44620cf1dd8SToby Isaac ierr = (*sp->ops->duplicate)(sp, spNew);CHKERRQ(ierr); 44720cf1dd8SToby Isaac PetscFunctionReturn(0); 44820cf1dd8SToby Isaac } 44920cf1dd8SToby Isaac 45020cf1dd8SToby Isaac /*@ 45120cf1dd8SToby Isaac PetscDualSpaceGetDM - Get the DM representing the reference cell 45220cf1dd8SToby Isaac 45320cf1dd8SToby Isaac Not collective 45420cf1dd8SToby Isaac 45520cf1dd8SToby Isaac Input Parameter: 45620cf1dd8SToby Isaac . sp - The PetscDualSpace 45720cf1dd8SToby Isaac 45820cf1dd8SToby Isaac Output Parameter: 45920cf1dd8SToby Isaac . dm - The reference cell 46020cf1dd8SToby Isaac 46120cf1dd8SToby Isaac Level: intermediate 46220cf1dd8SToby Isaac 46320cf1dd8SToby Isaac .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate() 46420cf1dd8SToby Isaac @*/ 46520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm) 46620cf1dd8SToby Isaac { 46720cf1dd8SToby Isaac PetscFunctionBegin; 46820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 46920cf1dd8SToby Isaac PetscValidPointer(dm, 2); 47020cf1dd8SToby Isaac *dm = sp->dm; 47120cf1dd8SToby Isaac PetscFunctionReturn(0); 47220cf1dd8SToby Isaac } 47320cf1dd8SToby Isaac 47420cf1dd8SToby Isaac /*@ 47520cf1dd8SToby Isaac PetscDualSpaceSetDM - Get the DM representing the reference cell 47620cf1dd8SToby Isaac 47720cf1dd8SToby Isaac Not collective 47820cf1dd8SToby Isaac 47920cf1dd8SToby Isaac Input Parameters: 48020cf1dd8SToby Isaac + sp - The PetscDualSpace 48120cf1dd8SToby Isaac - dm - The reference cell 48220cf1dd8SToby Isaac 48320cf1dd8SToby Isaac Level: intermediate 48420cf1dd8SToby Isaac 48520cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate() 48620cf1dd8SToby Isaac @*/ 48720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm) 48820cf1dd8SToby Isaac { 48920cf1dd8SToby Isaac PetscErrorCode ierr; 49020cf1dd8SToby Isaac 49120cf1dd8SToby Isaac PetscFunctionBegin; 49220cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 49320cf1dd8SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 49420cf1dd8SToby Isaac ierr = DMDestroy(&sp->dm);CHKERRQ(ierr); 49520cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 49620cf1dd8SToby Isaac sp->dm = dm; 49720cf1dd8SToby Isaac PetscFunctionReturn(0); 49820cf1dd8SToby Isaac } 49920cf1dd8SToby Isaac 50020cf1dd8SToby Isaac /*@ 50120cf1dd8SToby Isaac PetscDualSpaceGetOrder - Get the order of the dual space 50220cf1dd8SToby Isaac 50320cf1dd8SToby Isaac Not collective 50420cf1dd8SToby Isaac 50520cf1dd8SToby Isaac Input Parameter: 50620cf1dd8SToby Isaac . sp - The PetscDualSpace 50720cf1dd8SToby Isaac 50820cf1dd8SToby Isaac Output Parameter: 50920cf1dd8SToby Isaac . order - The order 51020cf1dd8SToby Isaac 51120cf1dd8SToby Isaac Level: intermediate 51220cf1dd8SToby Isaac 51320cf1dd8SToby Isaac .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate() 51420cf1dd8SToby Isaac @*/ 51520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order) 51620cf1dd8SToby Isaac { 51720cf1dd8SToby Isaac PetscFunctionBegin; 51820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 51920cf1dd8SToby Isaac PetscValidPointer(order, 2); 52020cf1dd8SToby Isaac *order = sp->order; 52120cf1dd8SToby Isaac PetscFunctionReturn(0); 52220cf1dd8SToby Isaac } 52320cf1dd8SToby Isaac 52420cf1dd8SToby Isaac /*@ 52520cf1dd8SToby Isaac PetscDualSpaceSetOrder - Set the order of the dual space 52620cf1dd8SToby Isaac 52720cf1dd8SToby Isaac Not collective 52820cf1dd8SToby Isaac 52920cf1dd8SToby Isaac Input Parameters: 53020cf1dd8SToby Isaac + sp - The PetscDualSpace 53120cf1dd8SToby Isaac - order - The order 53220cf1dd8SToby Isaac 53320cf1dd8SToby Isaac Level: intermediate 53420cf1dd8SToby Isaac 53520cf1dd8SToby Isaac .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate() 53620cf1dd8SToby Isaac @*/ 53720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order) 53820cf1dd8SToby Isaac { 53920cf1dd8SToby Isaac PetscFunctionBegin; 54020cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 54120cf1dd8SToby Isaac sp->order = order; 54220cf1dd8SToby Isaac PetscFunctionReturn(0); 54320cf1dd8SToby Isaac } 54420cf1dd8SToby Isaac 54520cf1dd8SToby Isaac /*@ 54620cf1dd8SToby Isaac PetscDualSpaceGetNumComponents - Return the number of components for this space 54720cf1dd8SToby Isaac 54820cf1dd8SToby Isaac Input Parameter: 54920cf1dd8SToby Isaac . sp - The PetscDualSpace 55020cf1dd8SToby Isaac 55120cf1dd8SToby Isaac Output Parameter: 55220cf1dd8SToby Isaac . Nc - The number of components 55320cf1dd8SToby Isaac 55420cf1dd8SToby Isaac Note: A vector space, for example, will have d components, where d is the spatial dimension 55520cf1dd8SToby Isaac 55620cf1dd8SToby Isaac Level: intermediate 55720cf1dd8SToby Isaac 55820cf1dd8SToby Isaac .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace 55920cf1dd8SToby Isaac @*/ 56020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc) 56120cf1dd8SToby Isaac { 56220cf1dd8SToby Isaac PetscFunctionBegin; 56320cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 56420cf1dd8SToby Isaac PetscValidPointer(Nc, 2); 56520cf1dd8SToby Isaac *Nc = sp->Nc; 56620cf1dd8SToby Isaac PetscFunctionReturn(0); 56720cf1dd8SToby Isaac } 56820cf1dd8SToby Isaac 56920cf1dd8SToby Isaac /*@ 57020cf1dd8SToby Isaac PetscDualSpaceSetNumComponents - Set the number of components for this space 57120cf1dd8SToby Isaac 57220cf1dd8SToby Isaac Input Parameters: 57320cf1dd8SToby Isaac + sp - The PetscDualSpace 57420cf1dd8SToby Isaac - order - The number of components 57520cf1dd8SToby Isaac 57620cf1dd8SToby Isaac Level: intermediate 57720cf1dd8SToby Isaac 57820cf1dd8SToby Isaac .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace 57920cf1dd8SToby Isaac @*/ 58020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc) 58120cf1dd8SToby Isaac { 58220cf1dd8SToby Isaac PetscFunctionBegin; 58320cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 58420cf1dd8SToby Isaac sp->Nc = Nc; 58520cf1dd8SToby Isaac PetscFunctionReturn(0); 58620cf1dd8SToby Isaac } 58720cf1dd8SToby Isaac 58820cf1dd8SToby Isaac /*@ 58920cf1dd8SToby Isaac PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space 59020cf1dd8SToby Isaac 59120cf1dd8SToby Isaac Not collective 59220cf1dd8SToby Isaac 59320cf1dd8SToby Isaac Input Parameters: 59420cf1dd8SToby Isaac + sp - The PetscDualSpace 59520cf1dd8SToby Isaac - i - The basis number 59620cf1dd8SToby Isaac 59720cf1dd8SToby Isaac Output Parameter: 59820cf1dd8SToby Isaac . functional - The basis functional 59920cf1dd8SToby Isaac 60020cf1dd8SToby Isaac Level: intermediate 60120cf1dd8SToby Isaac 60220cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate() 60320cf1dd8SToby Isaac @*/ 60420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional) 60520cf1dd8SToby Isaac { 60620cf1dd8SToby Isaac PetscInt dim; 60720cf1dd8SToby Isaac PetscErrorCode ierr; 60820cf1dd8SToby Isaac 60920cf1dd8SToby Isaac PetscFunctionBegin; 61020cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 61120cf1dd8SToby Isaac PetscValidPointer(functional, 3); 61220cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr); 61320cf1dd8SToby Isaac if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim); 61420cf1dd8SToby Isaac *functional = sp->functional[i]; 61520cf1dd8SToby Isaac PetscFunctionReturn(0); 61620cf1dd8SToby Isaac } 61720cf1dd8SToby Isaac 61820cf1dd8SToby Isaac /*@ 61920cf1dd8SToby Isaac PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals 62020cf1dd8SToby Isaac 62120cf1dd8SToby Isaac Not collective 62220cf1dd8SToby Isaac 62320cf1dd8SToby Isaac Input Parameter: 62420cf1dd8SToby Isaac . sp - The PetscDualSpace 62520cf1dd8SToby Isaac 62620cf1dd8SToby Isaac Output Parameter: 62720cf1dd8SToby Isaac . dim - The dimension 62820cf1dd8SToby Isaac 62920cf1dd8SToby Isaac Level: intermediate 63020cf1dd8SToby Isaac 63120cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 63220cf1dd8SToby Isaac @*/ 63320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim) 63420cf1dd8SToby Isaac { 63520cf1dd8SToby Isaac PetscErrorCode ierr; 63620cf1dd8SToby Isaac 63720cf1dd8SToby Isaac PetscFunctionBegin; 63820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 63920cf1dd8SToby Isaac PetscValidPointer(dim, 2); 64020cf1dd8SToby Isaac *dim = 0; 64120cf1dd8SToby Isaac if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);} 64220cf1dd8SToby Isaac PetscFunctionReturn(0); 64320cf1dd8SToby Isaac } 64420cf1dd8SToby Isaac 64520cf1dd8SToby Isaac /*@C 64620cf1dd8SToby Isaac PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension 64720cf1dd8SToby Isaac 64820cf1dd8SToby Isaac Not collective 64920cf1dd8SToby Isaac 65020cf1dd8SToby Isaac Input Parameter: 65120cf1dd8SToby Isaac . sp - The PetscDualSpace 65220cf1dd8SToby Isaac 65320cf1dd8SToby Isaac Output Parameter: 65420cf1dd8SToby Isaac . numDof - An array of length dim+1 which holds the number of dofs for each dimension 65520cf1dd8SToby Isaac 65620cf1dd8SToby Isaac Level: intermediate 65720cf1dd8SToby Isaac 65820cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 65920cf1dd8SToby Isaac @*/ 66020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof) 66120cf1dd8SToby Isaac { 66220cf1dd8SToby Isaac PetscErrorCode ierr; 66320cf1dd8SToby Isaac 66420cf1dd8SToby Isaac PetscFunctionBegin; 66520cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 66620cf1dd8SToby Isaac PetscValidPointer(numDof, 2); 66720cf1dd8SToby Isaac ierr = (*sp->ops->getnumdof)(sp, numDof);CHKERRQ(ierr); 66820cf1dd8SToby Isaac if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation"); 66920cf1dd8SToby Isaac PetscFunctionReturn(0); 67020cf1dd8SToby Isaac } 67120cf1dd8SToby Isaac 672a4ce7ad1SMatthew G. Knepley /*@ 673a4ce7ad1SMatthew G. Knepley PetscDualSpaceCreateSection - Create a PetscSection over the reference cell with the layout from this space 674a4ce7ad1SMatthew G. Knepley 675a4ce7ad1SMatthew G. Knepley Collective on sp 676a4ce7ad1SMatthew G. Knepley 677a4ce7ad1SMatthew G. Knepley Input Parameters: 678*f0fc11ceSJed Brown . sp - The PetscDualSpace 679a4ce7ad1SMatthew G. Knepley 680a4ce7ad1SMatthew G. Knepley Output Parameter: 681a4ce7ad1SMatthew G. Knepley . section - The section 682a4ce7ad1SMatthew G. Knepley 683a4ce7ad1SMatthew G. Knepley Level: advanced 684a4ce7ad1SMatthew G. Knepley 685a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate(), DMPLEX 686a4ce7ad1SMatthew G. Knepley @*/ 68720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace sp, PetscSection *section) 68820cf1dd8SToby Isaac { 68920cf1dd8SToby Isaac DM dm; 69020cf1dd8SToby Isaac PetscInt pStart, pEnd, depth, h, offset; 69120cf1dd8SToby Isaac const PetscInt *numDof; 69220cf1dd8SToby Isaac PetscErrorCode ierr; 69320cf1dd8SToby Isaac 69420cf1dd8SToby Isaac PetscFunctionBegin; 69520cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 69620cf1dd8SToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 69720cf1dd8SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sp),section);CHKERRQ(ierr); 69820cf1dd8SToby Isaac ierr = PetscSectionSetChart(*section,pStart,pEnd);CHKERRQ(ierr); 69920cf1dd8SToby Isaac ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 70020cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumDof(sp,&numDof);CHKERRQ(ierr); 70120cf1dd8SToby Isaac for (h = 0; h <= depth; h++) { 70220cf1dd8SToby Isaac PetscInt hStart, hEnd, p, dof; 70320cf1dd8SToby Isaac 70420cf1dd8SToby Isaac ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr); 70520cf1dd8SToby Isaac dof = numDof[depth - h]; 70620cf1dd8SToby Isaac for (p = hStart; p < hEnd; p++) { 70720cf1dd8SToby Isaac ierr = PetscSectionSetDof(*section,p,dof);CHKERRQ(ierr); 70820cf1dd8SToby Isaac } 70920cf1dd8SToby Isaac } 71020cf1dd8SToby Isaac ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 71120cf1dd8SToby Isaac for (h = 0, offset = 0; h <= depth; h++) { 71220cf1dd8SToby Isaac PetscInt hStart, hEnd, p, dof; 71320cf1dd8SToby Isaac 71420cf1dd8SToby Isaac ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr); 71520cf1dd8SToby Isaac dof = numDof[depth - h]; 71620cf1dd8SToby Isaac for (p = hStart; p < hEnd; p++) { 71720cf1dd8SToby Isaac ierr = PetscSectionGetDof(*section,p,&dof);CHKERRQ(ierr); 71820cf1dd8SToby Isaac ierr = PetscSectionSetOffset(*section,p,offset);CHKERRQ(ierr); 71920cf1dd8SToby Isaac offset += dof; 72020cf1dd8SToby Isaac } 72120cf1dd8SToby Isaac } 72220cf1dd8SToby Isaac PetscFunctionReturn(0); 72320cf1dd8SToby Isaac } 72420cf1dd8SToby Isaac 72520cf1dd8SToby Isaac /*@ 72620cf1dd8SToby Isaac PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell 72720cf1dd8SToby Isaac 728d083f849SBarry Smith Collective on sp 72920cf1dd8SToby Isaac 73020cf1dd8SToby Isaac Input Parameters: 73120cf1dd8SToby Isaac + sp - The PetscDualSpace 73220cf1dd8SToby Isaac . dim - The spatial dimension 73320cf1dd8SToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell 73420cf1dd8SToby Isaac 73520cf1dd8SToby Isaac Output Parameter: 73620cf1dd8SToby Isaac . refdm - The reference cell 73720cf1dd8SToby Isaac 738a4ce7ad1SMatthew G. Knepley Level: intermediate 73920cf1dd8SToby Isaac 74020cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), DMPLEX 74120cf1dd8SToby Isaac @*/ 74220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm) 74320cf1dd8SToby Isaac { 74420cf1dd8SToby Isaac PetscErrorCode ierr; 74520cf1dd8SToby Isaac 74620cf1dd8SToby Isaac PetscFunctionBeginUser; 74720cf1dd8SToby Isaac ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr); 74820cf1dd8SToby Isaac PetscFunctionReturn(0); 74920cf1dd8SToby Isaac } 75020cf1dd8SToby Isaac 75120cf1dd8SToby Isaac /*@C 75220cf1dd8SToby Isaac PetscDualSpaceApply - Apply a functional from the dual space basis to an input function 75320cf1dd8SToby Isaac 75420cf1dd8SToby Isaac Input Parameters: 75520cf1dd8SToby Isaac + sp - The PetscDualSpace object 75620cf1dd8SToby Isaac . f - The basis functional index 75720cf1dd8SToby Isaac . time - The time 75820cf1dd8SToby 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) 75920cf1dd8SToby Isaac . numComp - The number of components for the function 76020cf1dd8SToby Isaac . func - The input function 76120cf1dd8SToby Isaac - ctx - A context for the function 76220cf1dd8SToby Isaac 76320cf1dd8SToby Isaac Output Parameter: 76420cf1dd8SToby Isaac . value - numComp output values 76520cf1dd8SToby Isaac 76620cf1dd8SToby Isaac Note: The calling sequence for the callback func is given by: 76720cf1dd8SToby Isaac 76820cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[], 76920cf1dd8SToby Isaac $ PetscInt numComponents, PetscScalar values[], void *ctx) 77020cf1dd8SToby Isaac 771a4ce7ad1SMatthew G. Knepley Level: beginner 77220cf1dd8SToby Isaac 77320cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate() 77420cf1dd8SToby Isaac @*/ 77520cf1dd8SToby 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) 77620cf1dd8SToby Isaac { 77720cf1dd8SToby Isaac PetscErrorCode ierr; 77820cf1dd8SToby Isaac 77920cf1dd8SToby Isaac PetscFunctionBegin; 78020cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 78120cf1dd8SToby Isaac PetscValidPointer(cgeom, 4); 78220cf1dd8SToby Isaac PetscValidPointer(value, 8); 78320cf1dd8SToby Isaac ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr); 78420cf1dd8SToby Isaac PetscFunctionReturn(0); 78520cf1dd8SToby Isaac } 78620cf1dd8SToby Isaac 78720cf1dd8SToby Isaac /*@C 78820cf1dd8SToby Isaac PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints() 78920cf1dd8SToby Isaac 79020cf1dd8SToby Isaac Input Parameters: 79120cf1dd8SToby Isaac + sp - The PetscDualSpace object 79220cf1dd8SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints() 79320cf1dd8SToby Isaac 79420cf1dd8SToby Isaac Output Parameter: 79520cf1dd8SToby Isaac . spValue - The values of all dual space functionals 79620cf1dd8SToby Isaac 797a4ce7ad1SMatthew G. Knepley Level: beginner 79820cf1dd8SToby Isaac 79920cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate() 80020cf1dd8SToby Isaac @*/ 80120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 80220cf1dd8SToby Isaac { 80320cf1dd8SToby Isaac PetscErrorCode ierr; 80420cf1dd8SToby Isaac 80520cf1dd8SToby Isaac PetscFunctionBegin; 80620cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 80720cf1dd8SToby Isaac ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr); 80820cf1dd8SToby Isaac PetscFunctionReturn(0); 80920cf1dd8SToby Isaac } 81020cf1dd8SToby Isaac 81120cf1dd8SToby Isaac /*@C 81220cf1dd8SToby Isaac PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional. 81320cf1dd8SToby Isaac 81420cf1dd8SToby Isaac Input Parameters: 81520cf1dd8SToby Isaac + sp - The PetscDualSpace object 81620cf1dd8SToby Isaac . f - The basis functional index 81720cf1dd8SToby Isaac . time - The time 81820cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian) 81920cf1dd8SToby Isaac . Nc - The number of components for the function 82020cf1dd8SToby Isaac . func - The input function 82120cf1dd8SToby Isaac - ctx - A context for the function 82220cf1dd8SToby Isaac 82320cf1dd8SToby Isaac Output Parameter: 82420cf1dd8SToby Isaac . value - The output value 82520cf1dd8SToby Isaac 82620cf1dd8SToby Isaac Note: The calling sequence for the callback func is given by: 82720cf1dd8SToby Isaac 82820cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[], 82920cf1dd8SToby Isaac $ PetscInt numComponents, PetscScalar values[], void *ctx) 83020cf1dd8SToby Isaac 83120cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral 83220cf1dd8SToby Isaac 83320cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x) 83420cf1dd8SToby Isaac 83520cf1dd8SToby Isaac where both n and f have Nc components. 83620cf1dd8SToby Isaac 837a4ce7ad1SMatthew G. Knepley Level: beginner 83820cf1dd8SToby Isaac 83920cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate() 84020cf1dd8SToby Isaac @*/ 84120cf1dd8SToby 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) 84220cf1dd8SToby Isaac { 84320cf1dd8SToby Isaac DM dm; 84420cf1dd8SToby Isaac PetscQuadrature n; 84520cf1dd8SToby Isaac const PetscReal *points, *weights; 84620cf1dd8SToby Isaac PetscReal x[3]; 84720cf1dd8SToby Isaac PetscScalar *val; 84820cf1dd8SToby Isaac PetscInt dim, dE, qNc, c, Nq, q; 84920cf1dd8SToby Isaac PetscBool isAffine; 85020cf1dd8SToby Isaac PetscErrorCode ierr; 85120cf1dd8SToby Isaac 85220cf1dd8SToby Isaac PetscFunctionBegin; 85320cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 85420cf1dd8SToby Isaac PetscValidPointer(value, 5); 85520cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 85620cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 85720cf1dd8SToby Isaac ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 85820cf1dd8SToby Isaac if (dim != cgeom->dim) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim); 85920cf1dd8SToby Isaac if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 86020cf1dd8SToby Isaac ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 86120cf1dd8SToby Isaac *value = 0.0; 86220cf1dd8SToby Isaac isAffine = cgeom->isAffine; 86320cf1dd8SToby Isaac dE = cgeom->dimEmbed; 86420cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 86520cf1dd8SToby Isaac if (isAffine) { 86620cf1dd8SToby Isaac CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x); 86720cf1dd8SToby Isaac ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr); 86820cf1dd8SToby Isaac } else { 86920cf1dd8SToby Isaac ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr); 87020cf1dd8SToby Isaac } 87120cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 87220cf1dd8SToby Isaac *value += val[c]*weights[q*Nc+c]; 87320cf1dd8SToby Isaac } 87420cf1dd8SToby Isaac } 87520cf1dd8SToby Isaac ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 87620cf1dd8SToby Isaac PetscFunctionReturn(0); 87720cf1dd8SToby Isaac } 87820cf1dd8SToby Isaac 87920cf1dd8SToby Isaac /*@C 88020cf1dd8SToby Isaac PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints() 88120cf1dd8SToby Isaac 88220cf1dd8SToby Isaac Input Parameters: 88320cf1dd8SToby Isaac + sp - The PetscDualSpace object 88420cf1dd8SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints() 88520cf1dd8SToby Isaac 88620cf1dd8SToby Isaac Output Parameter: 88720cf1dd8SToby Isaac . spValue - The values of all dual space functionals 88820cf1dd8SToby Isaac 889a4ce7ad1SMatthew G. Knepley Level: beginner 89020cf1dd8SToby Isaac 89120cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate() 89220cf1dd8SToby Isaac @*/ 89320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 89420cf1dd8SToby Isaac { 89520cf1dd8SToby Isaac PetscQuadrature n; 89620cf1dd8SToby Isaac const PetscReal *points, *weights; 89720cf1dd8SToby Isaac PetscInt qNc, c, Nq, q, f, spdim, Nc; 89820cf1dd8SToby Isaac PetscInt offset; 89920cf1dd8SToby Isaac PetscErrorCode ierr; 90020cf1dd8SToby Isaac 90120cf1dd8SToby Isaac PetscFunctionBegin; 90220cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 90320cf1dd8SToby Isaac PetscValidScalarPointer(pointEval, 2); 90420cf1dd8SToby Isaac PetscValidScalarPointer(spValue, 5); 90520cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr); 90620cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 90720cf1dd8SToby Isaac for (f = 0, offset = 0; f < spdim; f++) { 90820cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 90920cf1dd8SToby Isaac ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 91020cf1dd8SToby Isaac if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 91120cf1dd8SToby Isaac spValue[f] = 0.0; 91220cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 91320cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 91420cf1dd8SToby Isaac spValue[f] += pointEval[offset++]*weights[q*Nc+c]; 91520cf1dd8SToby Isaac } 91620cf1dd8SToby Isaac } 91720cf1dd8SToby Isaac } 91820cf1dd8SToby Isaac PetscFunctionReturn(0); 91920cf1dd8SToby Isaac } 92020cf1dd8SToby Isaac 921a4ce7ad1SMatthew G. Knepley /*@ 922a4ce7ad1SMatthew G. Knepley PetscDualSpaceGetAllPoints - Get all quadrature points from this space 923a4ce7ad1SMatthew G. Knepley 924a4ce7ad1SMatthew G. Knepley Input Parameter: 925a4ce7ad1SMatthew G. Knepley . sp - The dualspace 926a4ce7ad1SMatthew G. Knepley 927a4ce7ad1SMatthew G. Knepley Output Parameter: 928a4ce7ad1SMatthew G. Knepley . allPoints - A PetscQuadrature object containing all evaluation points 929a4ce7ad1SMatthew G. Knepley 930a4ce7ad1SMatthew G. Knepley Level: advanced 931a4ce7ad1SMatthew G. Knepley 932a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate() 933a4ce7ad1SMatthew G. Knepley @*/ 93420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetAllPoints(PetscDualSpace sp, PetscQuadrature *allPoints) 93520cf1dd8SToby Isaac { 93620cf1dd8SToby Isaac PetscErrorCode ierr; 93720cf1dd8SToby Isaac 93820cf1dd8SToby Isaac PetscFunctionBegin; 93920cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 94020cf1dd8SToby Isaac PetscValidPointer(allPoints,2); 94120cf1dd8SToby Isaac if (!sp->allPoints && sp->ops->createallpoints) { 94220cf1dd8SToby Isaac ierr = (*sp->ops->createallpoints)(sp,&sp->allPoints);CHKERRQ(ierr); 94320cf1dd8SToby Isaac } 94420cf1dd8SToby Isaac *allPoints = sp->allPoints; 94520cf1dd8SToby Isaac PetscFunctionReturn(0); 94620cf1dd8SToby Isaac } 94720cf1dd8SToby Isaac 948a4ce7ad1SMatthew G. Knepley /*@ 949a4ce7ad1SMatthew G. Knepley PetscDualSpaceCreateAllPointsDefault - Create all evaluation points by examining functionals 950a4ce7ad1SMatthew G. Knepley 951a4ce7ad1SMatthew G. Knepley Input Parameter: 952a4ce7ad1SMatthew G. Knepley . sp - The dualspace 953a4ce7ad1SMatthew G. Knepley 954a4ce7ad1SMatthew G. Knepley Output Parameter: 955a4ce7ad1SMatthew G. Knepley . allPoints - A PetscQuadrature object containing all evaluation points 956a4ce7ad1SMatthew G. Knepley 957a4ce7ad1SMatthew G. Knepley Level: advanced 958a4ce7ad1SMatthew G. Knepley 959a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate() 960a4ce7ad1SMatthew G. Knepley @*/ 96120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateAllPointsDefault(PetscDualSpace sp, PetscQuadrature *allPoints) 96220cf1dd8SToby Isaac { 96320cf1dd8SToby Isaac PetscInt spdim; 96420cf1dd8SToby Isaac PetscInt numPoints, offset; 96520cf1dd8SToby Isaac PetscReal *points; 96620cf1dd8SToby Isaac PetscInt f, dim; 96720cf1dd8SToby Isaac PetscQuadrature q; 96820cf1dd8SToby Isaac PetscErrorCode ierr; 96920cf1dd8SToby Isaac 97020cf1dd8SToby Isaac PetscFunctionBegin; 97120cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr); 97220cf1dd8SToby Isaac if (!spdim) { 97320cf1dd8SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 97420cf1dd8SToby Isaac ierr = PetscQuadratureSetData(*allPoints,0,0,0,NULL,NULL);CHKERRQ(ierr); 97520cf1dd8SToby Isaac } 97620cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr); 97720cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr); 97820cf1dd8SToby Isaac for (f = 1; f < spdim; f++) { 97920cf1dd8SToby Isaac PetscInt Np; 98020cf1dd8SToby Isaac 98120cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr); 98220cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr); 98320cf1dd8SToby Isaac numPoints += Np; 98420cf1dd8SToby Isaac } 98520cf1dd8SToby Isaac ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 98620cf1dd8SToby Isaac for (f = 0, offset = 0; f < spdim; f++) { 98720cf1dd8SToby Isaac const PetscReal *p; 98820cf1dd8SToby Isaac PetscInt Np, i; 98920cf1dd8SToby Isaac 99020cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr); 99120cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,NULL);CHKERRQ(ierr); 99220cf1dd8SToby Isaac for (i = 0; i < Np * dim; i++) { 99320cf1dd8SToby Isaac points[offset + i] = p[i]; 99420cf1dd8SToby Isaac } 99520cf1dd8SToby Isaac offset += Np * dim; 99620cf1dd8SToby Isaac } 99720cf1dd8SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 99820cf1dd8SToby Isaac ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 99920cf1dd8SToby Isaac PetscFunctionReturn(0); 100020cf1dd8SToby Isaac } 100120cf1dd8SToby Isaac 100220cf1dd8SToby Isaac /*@C 100320cf1dd8SToby Isaac PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid. 100420cf1dd8SToby Isaac 100520cf1dd8SToby Isaac Input Parameters: 100620cf1dd8SToby Isaac + sp - The PetscDualSpace object 100720cf1dd8SToby Isaac . f - The basis functional index 100820cf1dd8SToby Isaac . time - The time 100920cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we currently just use the centroid 101020cf1dd8SToby Isaac . Nc - The number of components for the function 101120cf1dd8SToby Isaac . func - The input function 101220cf1dd8SToby Isaac - ctx - A context for the function 101320cf1dd8SToby Isaac 101420cf1dd8SToby Isaac Output Parameter: 101520cf1dd8SToby Isaac . value - The output value (scalar) 101620cf1dd8SToby Isaac 101720cf1dd8SToby Isaac Note: The calling sequence for the callback func is given by: 101820cf1dd8SToby Isaac 101920cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[], 102020cf1dd8SToby Isaac $ PetscInt numComponents, PetscScalar values[], void *ctx) 102120cf1dd8SToby Isaac 102220cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral 102320cf1dd8SToby Isaac 102420cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x) 102520cf1dd8SToby Isaac 102620cf1dd8SToby Isaac where both n and f have Nc components. 102720cf1dd8SToby Isaac 1028a4ce7ad1SMatthew G. Knepley Level: beginner 102920cf1dd8SToby Isaac 103020cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate() 103120cf1dd8SToby Isaac @*/ 103220cf1dd8SToby 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) 103320cf1dd8SToby Isaac { 103420cf1dd8SToby Isaac DM dm; 103520cf1dd8SToby Isaac PetscQuadrature n; 103620cf1dd8SToby Isaac const PetscReal *points, *weights; 103720cf1dd8SToby Isaac PetscScalar *val; 103820cf1dd8SToby Isaac PetscInt dimEmbed, qNc, c, Nq, q; 103920cf1dd8SToby Isaac PetscErrorCode ierr; 104020cf1dd8SToby Isaac 104120cf1dd8SToby Isaac PetscFunctionBegin; 104220cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 104320cf1dd8SToby Isaac PetscValidPointer(value, 5); 104420cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 104520cf1dd8SToby Isaac ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 104620cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 104720cf1dd8SToby Isaac ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 104820cf1dd8SToby Isaac if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 104920cf1dd8SToby Isaac ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 105020cf1dd8SToby Isaac *value = 0.; 105120cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 105220cf1dd8SToby Isaac ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr); 105320cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 105420cf1dd8SToby Isaac *value += val[c]*weights[q*Nc+c]; 105520cf1dd8SToby Isaac } 105620cf1dd8SToby Isaac } 105720cf1dd8SToby Isaac ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 105820cf1dd8SToby Isaac PetscFunctionReturn(0); 105920cf1dd8SToby Isaac } 106020cf1dd8SToby Isaac 106120cf1dd8SToby Isaac /*@ 106220cf1dd8SToby Isaac PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a 106320cf1dd8SToby Isaac given height. This assumes that the reference cell is symmetric over points of this height. 106420cf1dd8SToby Isaac 106520cf1dd8SToby Isaac If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and 106620cf1dd8SToby Isaac pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not 106720cf1dd8SToby Isaac support extracting subspaces, then NULL is returned. 106820cf1dd8SToby Isaac 106920cf1dd8SToby Isaac This does not increment the reference count on the returned dual space, and the user should not destroy it. 107020cf1dd8SToby Isaac 107120cf1dd8SToby Isaac Not collective 107220cf1dd8SToby Isaac 107320cf1dd8SToby Isaac Input Parameters: 107420cf1dd8SToby Isaac + sp - the PetscDualSpace object 107520cf1dd8SToby Isaac - height - the height of the mesh point for which the subspace is desired 107620cf1dd8SToby Isaac 107720cf1dd8SToby Isaac Output Parameter: 107820cf1dd8SToby Isaac . subsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the 107920cf1dd8SToby Isaac point, which will be of lesser dimension if height > 0. 108020cf1dd8SToby Isaac 108120cf1dd8SToby Isaac Level: advanced 108220cf1dd8SToby Isaac 108320cf1dd8SToby Isaac .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace 108420cf1dd8SToby Isaac @*/ 108520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp) 108620cf1dd8SToby Isaac { 108720cf1dd8SToby Isaac PetscErrorCode ierr; 108820cf1dd8SToby Isaac 108920cf1dd8SToby Isaac PetscFunctionBegin; 109020cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 109120cf1dd8SToby Isaac PetscValidPointer(subsp, 3); 109220cf1dd8SToby Isaac *subsp = NULL; 1093aa545514SMatthew G. Knepley if (sp->ops->getheightsubspace) {ierr = (*sp->ops->getheightsubspace)(sp, height, subsp);CHKERRQ(ierr);} 109420cf1dd8SToby Isaac PetscFunctionReturn(0); 109520cf1dd8SToby Isaac } 109620cf1dd8SToby Isaac 109720cf1dd8SToby Isaac /*@ 109820cf1dd8SToby Isaac PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point. 109920cf1dd8SToby Isaac 110020cf1dd8SToby Isaac If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not 110120cf1dd8SToby Isaac defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting 110220cf1dd8SToby Isaac subspaces, then NULL is returned. 110320cf1dd8SToby Isaac 110420cf1dd8SToby Isaac This does not increment the reference count on the returned dual space, and the user should not destroy it. 110520cf1dd8SToby Isaac 110620cf1dd8SToby Isaac Not collective 110720cf1dd8SToby Isaac 110820cf1dd8SToby Isaac Input Parameters: 110920cf1dd8SToby Isaac + sp - the PetscDualSpace object 111020cf1dd8SToby Isaac - point - the point (in the dual space's DM) for which the subspace is desired 111120cf1dd8SToby Isaac 111220cf1dd8SToby Isaac Output Parameters: 111320cf1dd8SToby Isaac bdsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the 111420cf1dd8SToby Isaac point, which will be of lesser dimension if height > 0. 111520cf1dd8SToby Isaac 111620cf1dd8SToby Isaac Level: advanced 111720cf1dd8SToby Isaac 111820cf1dd8SToby Isaac .seealso: PetscDualSpace 111920cf1dd8SToby Isaac @*/ 112020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp) 112120cf1dd8SToby Isaac { 112220cf1dd8SToby Isaac PetscErrorCode ierr; 112320cf1dd8SToby Isaac 112420cf1dd8SToby Isaac PetscFunctionBegin; 112520cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 112620cf1dd8SToby Isaac PetscValidPointer(bdsp,2); 112720cf1dd8SToby Isaac *bdsp = NULL; 112820cf1dd8SToby Isaac if (sp->ops->getpointsubspace) { 112920cf1dd8SToby Isaac ierr = (*sp->ops->getpointsubspace)(sp,point,bdsp);CHKERRQ(ierr); 113020cf1dd8SToby Isaac } else if (sp->ops->getheightsubspace) { 113120cf1dd8SToby Isaac DM dm; 113220cf1dd8SToby Isaac DMLabel label; 113320cf1dd8SToby Isaac PetscInt dim, depth, height; 113420cf1dd8SToby Isaac 113520cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 113620cf1dd8SToby Isaac ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr); 113720cf1dd8SToby Isaac ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr); 113820cf1dd8SToby Isaac ierr = DMLabelGetValue(label,point,&depth);CHKERRQ(ierr); 113920cf1dd8SToby Isaac height = dim - depth; 114020cf1dd8SToby Isaac ierr = (*sp->ops->getheightsubspace)(sp,height,bdsp);CHKERRQ(ierr); 114120cf1dd8SToby Isaac } 114220cf1dd8SToby Isaac PetscFunctionReturn(0); 114320cf1dd8SToby Isaac } 114420cf1dd8SToby Isaac 11456f905325SMatthew G. Knepley /*@C 11466f905325SMatthew G. Knepley PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis 11476f905325SMatthew G. Knepley 11486f905325SMatthew G. Knepley Not collective 11496f905325SMatthew G. Knepley 11506f905325SMatthew G. Knepley Input Parameter: 11516f905325SMatthew G. Knepley . sp - the PetscDualSpace object 11526f905325SMatthew G. Knepley 11536f905325SMatthew G. Knepley Output Parameters: 11546f905325SMatthew G. Knepley + perms - Permutations of the local degrees of freedom, parameterized by the point orientation 11556f905325SMatthew G. Knepley - flips - Sign reversal of the local degrees of freedom, parameterized by the point orientation 11566f905325SMatthew G. Knepley 11576f905325SMatthew G. Knepley Note: The permutation and flip arrays are organized in the following way 11586f905325SMatthew G. Knepley $ perms[p][ornt][dof # on point] = new local dof # 11596f905325SMatthew G. Knepley $ flips[p][ornt][dof # on point] = reversal or not 11606f905325SMatthew G. Knepley 11616f905325SMatthew G. Knepley Level: developer 11626f905325SMatthew G. Knepley 11636f905325SMatthew G. Knepley .seealso: PetscDualSpaceSetSymmetries() 11646f905325SMatthew G. Knepley @*/ 11656f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) 11666f905325SMatthew G. Knepley { 11676f905325SMatthew G. Knepley PetscErrorCode ierr; 11686f905325SMatthew G. Knepley 11696f905325SMatthew G. Knepley PetscFunctionBegin; 11706f905325SMatthew G. Knepley PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1); 11716f905325SMatthew G. Knepley if (perms) {PetscValidPointer(perms,2); *perms = NULL;} 11726f905325SMatthew G. Knepley if (flips) {PetscValidPointer(flips,3); *flips = NULL;} 11736f905325SMatthew G. Knepley if (sp->ops->getsymmetries) {ierr = (sp->ops->getsymmetries)(sp,perms,flips);CHKERRQ(ierr);} 11746f905325SMatthew G. Knepley PetscFunctionReturn(0); 11756f905325SMatthew G. Knepley } 11764bee2e38SMatthew G. Knepley 11774bee2e38SMatthew G. Knepley /*@ 11784bee2e38SMatthew G. Knepley PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space 11794bee2e38SMatthew G. Knepley 11804bee2e38SMatthew G. Knepley Input Parameter: 11814bee2e38SMatthew G. Knepley . dsp - The PetscDualSpace 11824bee2e38SMatthew G. Knepley 11834bee2e38SMatthew G. Knepley Output Parameter: 11844bee2e38SMatthew G. Knepley . k - The simplex dimension 11854bee2e38SMatthew G. Knepley 1186a4ce7ad1SMatthew G. Knepley Level: developer 11874bee2e38SMatthew G. Knepley 11884bee2e38SMatthew G. Knepley Note: Currently supported values are 11894bee2e38SMatthew G. Knepley $ 0: These are H_1 methods that only transform coordinates 11904bee2e38SMatthew G. Knepley $ 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM) 11914bee2e38SMatthew G. Knepley $ 2: These are the same as 1 11924bee2e38SMatthew G. Knepley $ 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM) 11934bee2e38SMatthew G. Knepley 11944bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType 11954bee2e38SMatthew G. Knepley @*/ 11964bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k) 11974bee2e38SMatthew G. Knepley { 11984bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 11994bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 12004bee2e38SMatthew G. Knepley PetscValidPointer(k, 2); 12014bee2e38SMatthew G. Knepley *k = dsp->k; 12024bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 12034bee2e38SMatthew G. Knepley } 12044bee2e38SMatthew G. Knepley 12054bee2e38SMatthew G. Knepley /*@C 12064bee2e38SMatthew G. Knepley PetscDualSpaceTransform - Transform the function values 12074bee2e38SMatthew G. Knepley 12084bee2e38SMatthew G. Knepley Input Parameters: 12094bee2e38SMatthew G. Knepley + dsp - The PetscDualSpace 12104bee2e38SMatthew G. Knepley . trans - The type of transform 12114bee2e38SMatthew G. Knepley . isInverse - Flag to invert the transform 12124bee2e38SMatthew G. Knepley . fegeom - The cell geometry 12134bee2e38SMatthew G. Knepley . Nv - The number of function samples 12144bee2e38SMatthew G. Knepley . Nc - The number of function components 12154bee2e38SMatthew G. Knepley - vals - The function values 12164bee2e38SMatthew G. Knepley 12174bee2e38SMatthew G. Knepley Output Parameter: 12184bee2e38SMatthew G. Knepley . vals - The transformed function values 12194bee2e38SMatthew G. Knepley 1220a4ce7ad1SMatthew G. Knepley Level: intermediate 12214bee2e38SMatthew G. Knepley 1222625e0966SMatthew G. Knepley .seealso: PetscDualSpaceTransformGradient(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType 12234bee2e38SMatthew G. Knepley @*/ 12244bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[]) 12254bee2e38SMatthew G. Knepley { 12264bee2e38SMatthew G. Knepley PetscInt dim, v, c; 12274bee2e38SMatthew G. Knepley 12284bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 12294bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 12304bee2e38SMatthew G. Knepley PetscValidPointer(fegeom, 4); 12314bee2e38SMatthew G. Knepley PetscValidPointer(vals, 7); 12322ae266adSMatthew G. Knepley dim = dsp->dm->dim; 12334bee2e38SMatthew G. Knepley /* Assume its a vector, otherwise assume its a bunch of scalars */ 12344bee2e38SMatthew G. Knepley if (Nc == 1 || Nc != dim) PetscFunctionReturn(0); 12354bee2e38SMatthew G. Knepley switch (trans) { 12364bee2e38SMatthew G. Knepley case IDENTITY_TRANSFORM: break; 12374bee2e38SMatthew G. Knepley case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */ 12384bee2e38SMatthew G. Knepley if (isInverse) { 12394bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 12404bee2e38SMatthew G. Knepley switch (dim) 12414bee2e38SMatthew G. Knepley { 12426142fa51SMatthew G. Knepley case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break; 12436142fa51SMatthew G. Knepley case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break; 12444bee2e38SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 12454bee2e38SMatthew G. Knepley } 12464bee2e38SMatthew G. Knepley } 12474bee2e38SMatthew G. Knepley } else { 12484bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 12494bee2e38SMatthew G. Knepley switch (dim) 12504bee2e38SMatthew G. Knepley { 12516142fa51SMatthew G. Knepley case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break; 12526142fa51SMatthew G. Knepley case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break; 12534bee2e38SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 12544bee2e38SMatthew G. Knepley } 12554bee2e38SMatthew G. Knepley } 12564bee2e38SMatthew G. Knepley } 12574bee2e38SMatthew G. Knepley break; 12584bee2e38SMatthew G. Knepley case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */ 12594bee2e38SMatthew G. Knepley if (isInverse) { 12604bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 12614bee2e38SMatthew G. Knepley switch (dim) 12624bee2e38SMatthew G. Knepley { 12636142fa51SMatthew G. Knepley case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break; 12646142fa51SMatthew G. Knepley case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break; 12654bee2e38SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 12664bee2e38SMatthew G. Knepley } 12674bee2e38SMatthew G. Knepley for (c = 0; c < Nc; ++c) vals[v*Nc+c] *= fegeom->detJ[0]; 12684bee2e38SMatthew G. Knepley } 12694bee2e38SMatthew G. Knepley } else { 12704bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 12714bee2e38SMatthew G. Knepley switch (dim) 12724bee2e38SMatthew G. Knepley { 12736142fa51SMatthew G. Knepley case 2: DMPlex_Mult2DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break; 12746142fa51SMatthew G. Knepley case 3: DMPlex_Mult3DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break; 12754bee2e38SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 12764bee2e38SMatthew G. Knepley } 12774bee2e38SMatthew G. Knepley for (c = 0; c < Nc; ++c) vals[v*Nc+c] /= fegeom->detJ[0]; 12784bee2e38SMatthew G. Knepley } 12794bee2e38SMatthew G. Knepley } 12804bee2e38SMatthew G. Knepley break; 12814bee2e38SMatthew G. Knepley } 12824bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 12834bee2e38SMatthew G. Knepley } 12844bee2e38SMatthew G. Knepley /*@C 12854bee2e38SMatthew G. Knepley PetscDualSpaceTransformGradient - Transform the function gradient values 12864bee2e38SMatthew G. Knepley 12874bee2e38SMatthew G. Knepley Input Parameters: 12884bee2e38SMatthew G. Knepley + dsp - The PetscDualSpace 12894bee2e38SMatthew G. Knepley . trans - The type of transform 12904bee2e38SMatthew G. Knepley . isInverse - Flag to invert the transform 12914bee2e38SMatthew G. Knepley . fegeom - The cell geometry 12924bee2e38SMatthew G. Knepley . Nv - The number of function gradient samples 12934bee2e38SMatthew G. Knepley . Nc - The number of function components 12944bee2e38SMatthew G. Knepley - vals - The function gradient values 12954bee2e38SMatthew G. Knepley 12964bee2e38SMatthew G. Knepley Output Parameter: 12974bee2e38SMatthew G. Knepley . vals - The transformed function values 12984bee2e38SMatthew G. Knepley 1299a4ce7ad1SMatthew G. Knepley Level: intermediate 13004bee2e38SMatthew G. Knepley 1301625e0966SMatthew G. Knepley .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType 13024bee2e38SMatthew G. Knepley @*/ 13034bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[]) 13044bee2e38SMatthew G. Knepley { 13054bee2e38SMatthew G. Knepley PetscInt dim, v, c, d; 13064bee2e38SMatthew G. Knepley 13074bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 13084bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 13094bee2e38SMatthew G. Knepley PetscValidPointer(fegeom, 4); 13104bee2e38SMatthew G. Knepley PetscValidPointer(vals, 7); 13112ae266adSMatthew G. Knepley dim = dsp->dm->dim; 13124bee2e38SMatthew G. Knepley /* Transform gradient */ 13134bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 13144bee2e38SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 13154bee2e38SMatthew G. Knepley switch (dim) 13164bee2e38SMatthew G. Knepley { 1317100a78e1SStefano Zampini case 1: vals[(v*Nc+c)*dim] *= fegeom->invJ[0];break; 13186142fa51SMatthew G. Knepley case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break; 13196142fa51SMatthew G. Knepley case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break; 13204bee2e38SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 13214bee2e38SMatthew G. Knepley } 13224bee2e38SMatthew G. Knepley } 13234bee2e38SMatthew G. Knepley } 13244bee2e38SMatthew G. Knepley /* Assume its a vector, otherwise assume its a bunch of scalars */ 13254bee2e38SMatthew G. Knepley if (Nc == 1 || Nc != dim) PetscFunctionReturn(0); 13264bee2e38SMatthew G. Knepley switch (trans) { 13274bee2e38SMatthew G. Knepley case IDENTITY_TRANSFORM: break; 13284bee2e38SMatthew G. Knepley case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */ 13294bee2e38SMatthew G. Knepley if (isInverse) { 13304bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 13314bee2e38SMatthew G. Knepley for (d = 0; d < dim; ++d) { 13324bee2e38SMatthew G. Knepley switch (dim) 13334bee2e38SMatthew G. Knepley { 13346142fa51SMatthew G. Knepley case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 13356142fa51SMatthew G. Knepley case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 13364bee2e38SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 13374bee2e38SMatthew G. Knepley } 13384bee2e38SMatthew G. Knepley } 13394bee2e38SMatthew G. Knepley } 13404bee2e38SMatthew G. Knepley } else { 13414bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 13424bee2e38SMatthew G. Knepley for (d = 0; d < dim; ++d) { 13434bee2e38SMatthew G. Knepley switch (dim) 13444bee2e38SMatthew G. Knepley { 13456142fa51SMatthew G. Knepley case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 13466142fa51SMatthew G. Knepley case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 13474bee2e38SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 13484bee2e38SMatthew G. Knepley } 13494bee2e38SMatthew G. Knepley } 13504bee2e38SMatthew G. Knepley } 13514bee2e38SMatthew G. Knepley } 13524bee2e38SMatthew G. Knepley break; 13534bee2e38SMatthew G. Knepley case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */ 13544bee2e38SMatthew G. Knepley if (isInverse) { 13554bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 13564bee2e38SMatthew G. Knepley for (d = 0; d < dim; ++d) { 13574bee2e38SMatthew G. Knepley switch (dim) 13584bee2e38SMatthew G. Knepley { 13596142fa51SMatthew G. Knepley case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 13606142fa51SMatthew G. Knepley case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 13614bee2e38SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 13624bee2e38SMatthew G. Knepley } 13634bee2e38SMatthew G. Knepley for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] *= fegeom->detJ[0]; 13644bee2e38SMatthew G. Knepley } 13654bee2e38SMatthew G. Knepley } 13664bee2e38SMatthew G. Knepley } else { 13674bee2e38SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 13684bee2e38SMatthew G. Knepley for (d = 0; d < dim; ++d) { 13694bee2e38SMatthew G. Knepley switch (dim) 13704bee2e38SMatthew G. Knepley { 13716142fa51SMatthew G. Knepley case 2: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 13726142fa51SMatthew G. Knepley case 3: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 13734bee2e38SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 13744bee2e38SMatthew G. Knepley } 13754bee2e38SMatthew G. Knepley for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] /= fegeom->detJ[0]; 13764bee2e38SMatthew G. Knepley } 13774bee2e38SMatthew G. Knepley } 13784bee2e38SMatthew G. Knepley } 13794bee2e38SMatthew G. Knepley break; 13804bee2e38SMatthew G. Knepley } 13814bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 13824bee2e38SMatthew G. Knepley } 13834bee2e38SMatthew G. Knepley 13844bee2e38SMatthew G. Knepley /*@C 13854bee2e38SMatthew 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. 13864bee2e38SMatthew G. Knepley 13874bee2e38SMatthew G. Knepley Input Parameters: 13884bee2e38SMatthew G. Knepley + dsp - The PetscDualSpace 13894bee2e38SMatthew G. Knepley . fegeom - The geometry for this cell 13904bee2e38SMatthew G. Knepley . Nq - The number of function samples 13914bee2e38SMatthew G. Knepley . Nc - The number of function components 13924bee2e38SMatthew G. Knepley - pointEval - The function values 13934bee2e38SMatthew G. Knepley 13944bee2e38SMatthew G. Knepley Output Parameter: 13954bee2e38SMatthew G. Knepley . pointEval - The transformed function values 13964bee2e38SMatthew G. Knepley 13974bee2e38SMatthew G. Knepley Level: advanced 13984bee2e38SMatthew G. Knepley 13994bee2e38SMatthew 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. 14004bee2e38SMatthew G. Knepley 14014bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 14024bee2e38SMatthew G. Knepley @*/ 14034bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 14044bee2e38SMatthew G. Knepley { 14054bee2e38SMatthew G. Knepley PetscDualSpaceTransformType trans; 14064bee2e38SMatthew G. Knepley PetscErrorCode ierr; 14074bee2e38SMatthew G. Knepley 14084bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 14094bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 14104bee2e38SMatthew G. Knepley PetscValidPointer(fegeom, 2); 14114bee2e38SMatthew G. Knepley PetscValidPointer(pointEval, 5); 14124bee2e38SMatthew G. Knepley /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 14134bee2e38SMatthew G. Knepley This determines their transformation properties. */ 14142ae266adSMatthew G. Knepley switch (dsp->k) 14154bee2e38SMatthew G. Knepley { 14164bee2e38SMatthew G. Knepley case 0: /* H^1 point evaluations */ 14174bee2e38SMatthew G. Knepley trans = IDENTITY_TRANSFORM;break; 14184bee2e38SMatthew G. Knepley case 1: /* Hcurl preserves tangential edge traces */ 14194bee2e38SMatthew G. Knepley case 2: 14204bee2e38SMatthew G. Knepley trans = COVARIANT_PIOLA_TRANSFORM;break; 14214bee2e38SMatthew G. Knepley case 3: /* Hdiv preserve normal traces */ 14224bee2e38SMatthew G. Knepley trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 14232ae266adSMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k); 14244bee2e38SMatthew G. Knepley } 14254bee2e38SMatthew G. Knepley ierr = PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr); 14264bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 14274bee2e38SMatthew G. Knepley } 14284bee2e38SMatthew G. Knepley 14294bee2e38SMatthew G. Knepley /*@C 14304bee2e38SMatthew 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. 14314bee2e38SMatthew G. Knepley 14324bee2e38SMatthew G. Knepley Input Parameters: 14334bee2e38SMatthew G. Knepley + dsp - The PetscDualSpace 14344bee2e38SMatthew G. Knepley . fegeom - The geometry for this cell 14354bee2e38SMatthew G. Knepley . Nq - The number of function samples 14364bee2e38SMatthew G. Knepley . Nc - The number of function components 14374bee2e38SMatthew G. Knepley - pointEval - The function values 14384bee2e38SMatthew G. Knepley 14394bee2e38SMatthew G. Knepley Output Parameter: 14404bee2e38SMatthew G. Knepley . pointEval - The transformed function values 14414bee2e38SMatthew G. Knepley 14424bee2e38SMatthew G. Knepley Level: advanced 14434bee2e38SMatthew G. Knepley 14444bee2e38SMatthew 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. 14454bee2e38SMatthew G. Knepley 14464bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 14474bee2e38SMatthew G. Knepley @*/ 14484bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 14494bee2e38SMatthew G. Knepley { 14504bee2e38SMatthew G. Knepley PetscDualSpaceTransformType trans; 14514bee2e38SMatthew G. Knepley PetscErrorCode ierr; 14524bee2e38SMatthew G. Knepley 14534bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 14544bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 14554bee2e38SMatthew G. Knepley PetscValidPointer(fegeom, 2); 14564bee2e38SMatthew G. Knepley PetscValidPointer(pointEval, 5); 14574bee2e38SMatthew G. Knepley /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 14584bee2e38SMatthew G. Knepley This determines their transformation properties. */ 14592ae266adSMatthew G. Knepley switch (dsp->k) 14604bee2e38SMatthew G. Knepley { 14614bee2e38SMatthew G. Knepley case 0: /* H^1 point evaluations */ 14624bee2e38SMatthew G. Knepley trans = IDENTITY_TRANSFORM;break; 14634bee2e38SMatthew G. Knepley case 1: /* Hcurl preserves tangential edge traces */ 14644bee2e38SMatthew G. Knepley case 2: 14654bee2e38SMatthew G. Knepley trans = COVARIANT_PIOLA_TRANSFORM;break; 14664bee2e38SMatthew G. Knepley case 3: /* Hdiv preserve normal traces */ 14674bee2e38SMatthew G. Knepley trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 14682ae266adSMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k); 14694bee2e38SMatthew G. Knepley } 14704bee2e38SMatthew G. Knepley ierr = PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr); 14714bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 14724bee2e38SMatthew G. Knepley } 14734bee2e38SMatthew G. Knepley 14744bee2e38SMatthew G. Knepley /*@C 14754bee2e38SMatthew 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. 14764bee2e38SMatthew G. Knepley 14774bee2e38SMatthew G. Knepley Input Parameters: 14784bee2e38SMatthew G. Knepley + dsp - The PetscDualSpace 14794bee2e38SMatthew G. Knepley . fegeom - The geometry for this cell 14804bee2e38SMatthew G. Knepley . Nq - The number of function gradient samples 14814bee2e38SMatthew G. Knepley . Nc - The number of function components 14824bee2e38SMatthew G. Knepley - pointEval - The function gradient values 14834bee2e38SMatthew G. Knepley 14844bee2e38SMatthew G. Knepley Output Parameter: 14854bee2e38SMatthew G. Knepley . pointEval - The transformed function gradient values 14864bee2e38SMatthew G. Knepley 14874bee2e38SMatthew G. Knepley Level: advanced 14884bee2e38SMatthew G. Knepley 14894bee2e38SMatthew 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. 14904bee2e38SMatthew G. Knepley 14914bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 1492dc0529c6SBarry Smith @*/ 1493dc0529c6SBarry Smith PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 14944bee2e38SMatthew G. Knepley { 14954bee2e38SMatthew G. Knepley PetscDualSpaceTransformType trans; 14964bee2e38SMatthew G. Knepley PetscErrorCode ierr; 14974bee2e38SMatthew G. Knepley 14984bee2e38SMatthew G. Knepley PetscFunctionBeginHot; 14994bee2e38SMatthew G. Knepley PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 15004bee2e38SMatthew G. Knepley PetscValidPointer(fegeom, 2); 15014bee2e38SMatthew G. Knepley PetscValidPointer(pointEval, 5); 15024bee2e38SMatthew G. Knepley /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 15034bee2e38SMatthew G. Knepley This determines their transformation properties. */ 15042ae266adSMatthew G. Knepley switch (dsp->k) 15054bee2e38SMatthew G. Knepley { 15064bee2e38SMatthew G. Knepley case 0: /* H^1 point evaluations */ 15074bee2e38SMatthew G. Knepley trans = IDENTITY_TRANSFORM;break; 15084bee2e38SMatthew G. Knepley case 1: /* Hcurl preserves tangential edge traces */ 15094bee2e38SMatthew G. Knepley case 2: 15104bee2e38SMatthew G. Knepley trans = COVARIANT_PIOLA_TRANSFORM;break; 15114bee2e38SMatthew G. Knepley case 3: /* Hdiv preserve normal traces */ 15124bee2e38SMatthew G. Knepley trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 15132ae266adSMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k); 15144bee2e38SMatthew G. Knepley } 15154bee2e38SMatthew G. Knepley ierr = PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr); 15164bee2e38SMatthew G. Knepley PetscFunctionReturn(0); 15174bee2e38SMatthew G. Knepley } 1518