1c0ce001eSMatthew G. Knepley #include <petsc/private/petscfvimpl.h> /*I "petscfv.h" I*/ 2*412e9a14SMatthew G. Knepley #include <petscdmplex.h> 3c0ce001eSMatthew G. Knepley #include <petscds.h> 4c0ce001eSMatthew G. Knepley 5c0ce001eSMatthew G. Knepley PetscClassId PETSCLIMITER_CLASSID = 0; 6c0ce001eSMatthew G. Knepley 7c0ce001eSMatthew G. Knepley PetscFunctionList PetscLimiterList = NULL; 8c0ce001eSMatthew G. Knepley PetscBool PetscLimiterRegisterAllCalled = PETSC_FALSE; 9c0ce001eSMatthew G. Knepley 10c0ce001eSMatthew G. Knepley PetscBool Limitercite = PETSC_FALSE; 11c0ce001eSMatthew G. Knepley const char LimiterCitation[] = "@article{BergerAftosmisMurman2005,\n" 12c0ce001eSMatthew G. Knepley " title = {Analysis of slope limiters on irregular grids},\n" 13c0ce001eSMatthew G. Knepley " journal = {AIAA paper},\n" 14c0ce001eSMatthew G. Knepley " author = {Marsha Berger and Michael J. Aftosmis and Scott M. Murman},\n" 15c0ce001eSMatthew G. Knepley " volume = {490},\n" 16c0ce001eSMatthew G. Knepley " year = {2005}\n}\n"; 17c0ce001eSMatthew G. Knepley 18c0ce001eSMatthew G. Knepley /*@C 19c0ce001eSMatthew G. Knepley PetscLimiterRegister - Adds a new PetscLimiter implementation 20c0ce001eSMatthew G. Knepley 21c0ce001eSMatthew G. Knepley Not Collective 22c0ce001eSMatthew G. Knepley 23c0ce001eSMatthew G. Knepley Input Parameters: 24c0ce001eSMatthew G. Knepley + name - The name of a new user-defined creation routine 25c0ce001eSMatthew G. Knepley - create_func - The creation routine itself 26c0ce001eSMatthew G. Knepley 27c0ce001eSMatthew G. Knepley Notes: 28c0ce001eSMatthew G. Knepley PetscLimiterRegister() may be called multiple times to add several user-defined PetscLimiters 29c0ce001eSMatthew G. Knepley 30c0ce001eSMatthew G. Knepley Sample usage: 31c0ce001eSMatthew G. Knepley .vb 32c0ce001eSMatthew G. Knepley PetscLimiterRegister("my_lim", MyPetscLimiterCreate); 33c0ce001eSMatthew G. Knepley .ve 34c0ce001eSMatthew G. Knepley 35c0ce001eSMatthew G. Knepley Then, your PetscLimiter type can be chosen with the procedural interface via 36c0ce001eSMatthew G. Knepley .vb 37c0ce001eSMatthew G. Knepley PetscLimiterCreate(MPI_Comm, PetscLimiter *); 38c0ce001eSMatthew G. Knepley PetscLimiterSetType(PetscLimiter, "my_lim"); 39c0ce001eSMatthew G. Knepley .ve 40c0ce001eSMatthew G. Knepley or at runtime via the option 41c0ce001eSMatthew G. Knepley .vb 42c0ce001eSMatthew G. Knepley -petsclimiter_type my_lim 43c0ce001eSMatthew G. Knepley .ve 44c0ce001eSMatthew G. Knepley 45c0ce001eSMatthew G. Knepley Level: advanced 46c0ce001eSMatthew G. Knepley 47c0ce001eSMatthew G. Knepley .seealso: PetscLimiterRegisterAll(), PetscLimiterRegisterDestroy() 48c0ce001eSMatthew G. Knepley 49c0ce001eSMatthew G. Knepley @*/ 50c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterRegister(const char sname[], PetscErrorCode (*function)(PetscLimiter)) 51c0ce001eSMatthew G. Knepley { 52c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 53c0ce001eSMatthew G. Knepley 54c0ce001eSMatthew G. Knepley PetscFunctionBegin; 55c0ce001eSMatthew G. Knepley ierr = PetscFunctionListAdd(&PetscLimiterList, sname, function);CHKERRQ(ierr); 56c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 57c0ce001eSMatthew G. Knepley } 58c0ce001eSMatthew G. Knepley 59c0ce001eSMatthew G. Knepley /*@C 60c0ce001eSMatthew G. Knepley PetscLimiterSetType - Builds a particular PetscLimiter 61c0ce001eSMatthew G. Knepley 62c0ce001eSMatthew G. Knepley Collective on lim 63c0ce001eSMatthew G. Knepley 64c0ce001eSMatthew G. Knepley Input Parameters: 65c0ce001eSMatthew G. Knepley + lim - The PetscLimiter object 66c0ce001eSMatthew G. Knepley - name - The kind of limiter 67c0ce001eSMatthew G. Knepley 68c0ce001eSMatthew G. Knepley Options Database Key: 69c0ce001eSMatthew G. Knepley . -petsclimiter_type <type> - Sets the PetscLimiter type; use -help for a list of available types 70c0ce001eSMatthew G. Knepley 71c0ce001eSMatthew G. Knepley Level: intermediate 72c0ce001eSMatthew G. Knepley 73c0ce001eSMatthew G. Knepley .seealso: PetscLimiterGetType(), PetscLimiterCreate() 74c0ce001eSMatthew G. Knepley @*/ 75c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterSetType(PetscLimiter lim, PetscLimiterType name) 76c0ce001eSMatthew G. Knepley { 77c0ce001eSMatthew G. Knepley PetscErrorCode (*r)(PetscLimiter); 78c0ce001eSMatthew G. Knepley PetscBool match; 79c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 80c0ce001eSMatthew G. Knepley 81c0ce001eSMatthew G. Knepley PetscFunctionBegin; 82c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 83c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) lim, name, &match);CHKERRQ(ierr); 84c0ce001eSMatthew G. Knepley if (match) PetscFunctionReturn(0); 85c0ce001eSMatthew G. Knepley 86c0ce001eSMatthew G. Knepley ierr = PetscLimiterRegisterAll();CHKERRQ(ierr); 87c0ce001eSMatthew G. Knepley ierr = PetscFunctionListFind(PetscLimiterList, name, &r);CHKERRQ(ierr); 88c0ce001eSMatthew G. Knepley if (!r) SETERRQ1(PetscObjectComm((PetscObject) lim), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscLimiter type: %s", name); 89c0ce001eSMatthew G. Knepley 90c0ce001eSMatthew G. Knepley if (lim->ops->destroy) { 91c0ce001eSMatthew G. Knepley ierr = (*lim->ops->destroy)(lim);CHKERRQ(ierr); 92c0ce001eSMatthew G. Knepley lim->ops->destroy = NULL; 93c0ce001eSMatthew G. Knepley } 94c0ce001eSMatthew G. Knepley ierr = (*r)(lim);CHKERRQ(ierr); 95c0ce001eSMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) lim, name);CHKERRQ(ierr); 96c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 97c0ce001eSMatthew G. Knepley } 98c0ce001eSMatthew G. Knepley 99c0ce001eSMatthew G. Knepley /*@C 100c0ce001eSMatthew G. Knepley PetscLimiterGetType - Gets the PetscLimiter type name (as a string) from the object. 101c0ce001eSMatthew G. Knepley 102c0ce001eSMatthew G. Knepley Not Collective 103c0ce001eSMatthew G. Knepley 104c0ce001eSMatthew G. Knepley Input Parameter: 105c0ce001eSMatthew G. Knepley . lim - The PetscLimiter 106c0ce001eSMatthew G. Knepley 107c0ce001eSMatthew G. Knepley Output Parameter: 108c0ce001eSMatthew G. Knepley . name - The PetscLimiter type name 109c0ce001eSMatthew G. Knepley 110c0ce001eSMatthew G. Knepley Level: intermediate 111c0ce001eSMatthew G. Knepley 112c0ce001eSMatthew G. Knepley .seealso: PetscLimiterSetType(), PetscLimiterCreate() 113c0ce001eSMatthew G. Knepley @*/ 114c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterGetType(PetscLimiter lim, PetscLimiterType *name) 115c0ce001eSMatthew G. Knepley { 116c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 117c0ce001eSMatthew G. Knepley 118c0ce001eSMatthew G. Knepley PetscFunctionBegin; 119c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 120c0ce001eSMatthew G. Knepley PetscValidPointer(name, 2); 121c0ce001eSMatthew G. Knepley ierr = PetscLimiterRegisterAll();CHKERRQ(ierr); 122c0ce001eSMatthew G. Knepley *name = ((PetscObject) lim)->type_name; 123c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 124c0ce001eSMatthew G. Knepley } 125c0ce001eSMatthew G. Knepley 126c0ce001eSMatthew G. Knepley /*@C 127fe2efc57SMark PetscLimiterViewFromOptions - View from Options 128fe2efc57SMark 129fe2efc57SMark Collective on PetscLimiter 130fe2efc57SMark 131fe2efc57SMark Input Parameters: 132fe2efc57SMark + A - the PetscLimiter object to view 133736c3998SJose E. Roman . obj - Optional object 134736c3998SJose E. Roman - name - command line option 135fe2efc57SMark 136fe2efc57SMark Level: intermediate 137fe2efc57SMark .seealso: PetscLimiter, PetscLimiterView, PetscObjectViewFromOptions(), PetscLimiterCreate() 138fe2efc57SMark @*/ 139fe2efc57SMark PetscErrorCode PetscLimiterViewFromOptions(PetscLimiter A,PetscObject obj,const char name[]) 140fe2efc57SMark { 141fe2efc57SMark PetscErrorCode ierr; 142fe2efc57SMark 143fe2efc57SMark PetscFunctionBegin; 144fe2efc57SMark PetscValidHeaderSpecific(A,PETSCLIMITER_CLASSID,1); 145fe2efc57SMark ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr); 146fe2efc57SMark PetscFunctionReturn(0); 147fe2efc57SMark } 148fe2efc57SMark 149fe2efc57SMark /*@C 150c0ce001eSMatthew G. Knepley PetscLimiterView - Views a PetscLimiter 151c0ce001eSMatthew G. Knepley 152c0ce001eSMatthew G. Knepley Collective on lim 153c0ce001eSMatthew G. Knepley 154c0ce001eSMatthew G. Knepley Input Parameter: 155c0ce001eSMatthew G. Knepley + lim - the PetscLimiter object to view 156c0ce001eSMatthew G. Knepley - v - the viewer 157c0ce001eSMatthew G. Knepley 15888f5f89eSMatthew G. Knepley Level: beginner 159c0ce001eSMatthew G. Knepley 160c0ce001eSMatthew G. Knepley .seealso: PetscLimiterDestroy() 161c0ce001eSMatthew G. Knepley @*/ 162c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView(PetscLimiter lim, PetscViewer v) 163c0ce001eSMatthew G. Knepley { 164c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 165c0ce001eSMatthew G. Knepley 166c0ce001eSMatthew G. Knepley PetscFunctionBegin; 167c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 168c0ce001eSMatthew G. Knepley if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) lim), &v);CHKERRQ(ierr);} 169c0ce001eSMatthew G. Knepley if (lim->ops->view) {ierr = (*lim->ops->view)(lim, v);CHKERRQ(ierr);} 170c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 171c0ce001eSMatthew G. Knepley } 172c0ce001eSMatthew G. Knepley 173c0ce001eSMatthew G. Knepley /*@ 174c0ce001eSMatthew G. Knepley PetscLimiterSetFromOptions - sets parameters in a PetscLimiter from the options database 175c0ce001eSMatthew G. Knepley 176c0ce001eSMatthew G. Knepley Collective on lim 177c0ce001eSMatthew G. Knepley 178c0ce001eSMatthew G. Knepley Input Parameter: 179c0ce001eSMatthew G. Knepley . lim - the PetscLimiter object to set options for 180c0ce001eSMatthew G. Knepley 18188f5f89eSMatthew G. Knepley Level: intermediate 182c0ce001eSMatthew G. Knepley 183c0ce001eSMatthew G. Knepley .seealso: PetscLimiterView() 184c0ce001eSMatthew G. Knepley @*/ 185c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterSetFromOptions(PetscLimiter lim) 186c0ce001eSMatthew G. Knepley { 187c0ce001eSMatthew G. Knepley const char *defaultType; 188c0ce001eSMatthew G. Knepley char name[256]; 189c0ce001eSMatthew G. Knepley PetscBool flg; 190c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 191c0ce001eSMatthew G. Knepley 192c0ce001eSMatthew G. Knepley PetscFunctionBegin; 193c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 194c0ce001eSMatthew G. Knepley if (!((PetscObject) lim)->type_name) defaultType = PETSCLIMITERSIN; 195c0ce001eSMatthew G. Knepley else defaultType = ((PetscObject) lim)->type_name; 196c0ce001eSMatthew G. Knepley ierr = PetscLimiterRegisterAll();CHKERRQ(ierr); 197c0ce001eSMatthew G. Knepley 198c0ce001eSMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) lim);CHKERRQ(ierr); 199c0ce001eSMatthew G. Knepley ierr = PetscOptionsFList("-petsclimiter_type", "Finite volume slope limiter", "PetscLimiterSetType", PetscLimiterList, defaultType, name, 256, &flg);CHKERRQ(ierr); 200c0ce001eSMatthew G. Knepley if (flg) { 201c0ce001eSMatthew G. Knepley ierr = PetscLimiterSetType(lim, name);CHKERRQ(ierr); 202c0ce001eSMatthew G. Knepley } else if (!((PetscObject) lim)->type_name) { 203c0ce001eSMatthew G. Knepley ierr = PetscLimiterSetType(lim, defaultType);CHKERRQ(ierr); 204c0ce001eSMatthew G. Knepley } 205c0ce001eSMatthew G. Knepley if (lim->ops->setfromoptions) {ierr = (*lim->ops->setfromoptions)(lim);CHKERRQ(ierr);} 206c0ce001eSMatthew G. Knepley /* process any options handlers added with PetscObjectAddOptionsHandler() */ 207c0ce001eSMatthew G. Knepley ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) lim);CHKERRQ(ierr); 208c0ce001eSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 209c0ce001eSMatthew G. Knepley ierr = PetscLimiterViewFromOptions(lim, NULL, "-petsclimiter_view");CHKERRQ(ierr); 210c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 211c0ce001eSMatthew G. Knepley } 212c0ce001eSMatthew G. Knepley 213c0ce001eSMatthew G. Knepley /*@C 214c0ce001eSMatthew G. Knepley PetscLimiterSetUp - Construct data structures for the PetscLimiter 215c0ce001eSMatthew G. Knepley 216c0ce001eSMatthew G. Knepley Collective on lim 217c0ce001eSMatthew G. Knepley 218c0ce001eSMatthew G. Knepley Input Parameter: 219c0ce001eSMatthew G. Knepley . lim - the PetscLimiter object to setup 220c0ce001eSMatthew G. Knepley 22188f5f89eSMatthew G. Knepley Level: intermediate 222c0ce001eSMatthew G. Knepley 223c0ce001eSMatthew G. Knepley .seealso: PetscLimiterView(), PetscLimiterDestroy() 224c0ce001eSMatthew G. Knepley @*/ 225c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterSetUp(PetscLimiter lim) 226c0ce001eSMatthew G. Knepley { 227c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 228c0ce001eSMatthew G. Knepley 229c0ce001eSMatthew G. Knepley PetscFunctionBegin; 230c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 231c0ce001eSMatthew G. Knepley if (lim->ops->setup) {ierr = (*lim->ops->setup)(lim);CHKERRQ(ierr);} 232c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 233c0ce001eSMatthew G. Knepley } 234c0ce001eSMatthew G. Knepley 235c0ce001eSMatthew G. Knepley /*@ 236c0ce001eSMatthew G. Knepley PetscLimiterDestroy - Destroys a PetscLimiter object 237c0ce001eSMatthew G. Knepley 238c0ce001eSMatthew G. Knepley Collective on lim 239c0ce001eSMatthew G. Knepley 240c0ce001eSMatthew G. Knepley Input Parameter: 241c0ce001eSMatthew G. Knepley . lim - the PetscLimiter object to destroy 242c0ce001eSMatthew G. Knepley 24388f5f89eSMatthew G. Knepley Level: beginner 244c0ce001eSMatthew G. Knepley 245c0ce001eSMatthew G. Knepley .seealso: PetscLimiterView() 246c0ce001eSMatthew G. Knepley @*/ 247c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterDestroy(PetscLimiter *lim) 248c0ce001eSMatthew G. Knepley { 249c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 250c0ce001eSMatthew G. Knepley 251c0ce001eSMatthew G. Knepley PetscFunctionBegin; 252c0ce001eSMatthew G. Knepley if (!*lim) PetscFunctionReturn(0); 253c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific((*lim), PETSCLIMITER_CLASSID, 1); 254c0ce001eSMatthew G. Knepley 255c0ce001eSMatthew G. Knepley if (--((PetscObject)(*lim))->refct > 0) {*lim = 0; PetscFunctionReturn(0);} 256c0ce001eSMatthew G. Knepley ((PetscObject) (*lim))->refct = 0; 257c0ce001eSMatthew G. Knepley 258c0ce001eSMatthew G. Knepley if ((*lim)->ops->destroy) {ierr = (*(*lim)->ops->destroy)(*lim);CHKERRQ(ierr);} 259c0ce001eSMatthew G. Knepley ierr = PetscHeaderDestroy(lim);CHKERRQ(ierr); 260c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 261c0ce001eSMatthew G. Knepley } 262c0ce001eSMatthew G. Knepley 263c0ce001eSMatthew G. Knepley /*@ 264c0ce001eSMatthew G. Knepley PetscLimiterCreate - Creates an empty PetscLimiter object. The type can then be set with PetscLimiterSetType(). 265c0ce001eSMatthew G. Knepley 266c0ce001eSMatthew G. Knepley Collective 267c0ce001eSMatthew G. Knepley 268c0ce001eSMatthew G. Knepley Input Parameter: 269c0ce001eSMatthew G. Knepley . comm - The communicator for the PetscLimiter object 270c0ce001eSMatthew G. Knepley 271c0ce001eSMatthew G. Knepley Output Parameter: 272c0ce001eSMatthew G. Knepley . lim - The PetscLimiter object 273c0ce001eSMatthew G. Knepley 274c0ce001eSMatthew G. Knepley Level: beginner 275c0ce001eSMatthew G. Knepley 276c0ce001eSMatthew G. Knepley .seealso: PetscLimiterSetType(), PETSCLIMITERSIN 277c0ce001eSMatthew G. Knepley @*/ 278c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterCreate(MPI_Comm comm, PetscLimiter *lim) 279c0ce001eSMatthew G. Knepley { 280c0ce001eSMatthew G. Knepley PetscLimiter l; 281c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 282c0ce001eSMatthew G. Knepley 283c0ce001eSMatthew G. Knepley PetscFunctionBegin; 284c0ce001eSMatthew G. Knepley PetscValidPointer(lim, 2); 285c0ce001eSMatthew G. Knepley ierr = PetscCitationsRegister(LimiterCitation,&Limitercite);CHKERRQ(ierr); 286c0ce001eSMatthew G. Knepley *lim = NULL; 287c0ce001eSMatthew G. Knepley ierr = PetscFVInitializePackage();CHKERRQ(ierr); 288c0ce001eSMatthew G. Knepley 289c0ce001eSMatthew G. Knepley ierr = PetscHeaderCreate(l, PETSCLIMITER_CLASSID, "PetscLimiter", "Finite Volume Slope Limiter", "PetscLimiter", comm, PetscLimiterDestroy, PetscLimiterView);CHKERRQ(ierr); 290c0ce001eSMatthew G. Knepley 291c0ce001eSMatthew G. Knepley *lim = l; 292c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 293c0ce001eSMatthew G. Knepley } 294c0ce001eSMatthew G. Knepley 29588f5f89eSMatthew G. Knepley /*@ 29688f5f89eSMatthew G. Knepley PetscLimiterLimit - Limit the flux 29788f5f89eSMatthew G. Knepley 29888f5f89eSMatthew G. Knepley Input Parameters: 29988f5f89eSMatthew G. Knepley + lim - The PetscLimiter 30088f5f89eSMatthew G. Knepley - flim - The input field 30188f5f89eSMatthew G. Knepley 30288f5f89eSMatthew G. Knepley Output Parameter: 30388f5f89eSMatthew G. Knepley . phi - The limited field 30488f5f89eSMatthew G. Knepley 30588f5f89eSMatthew G. Knepley Note: Limiters given in symmetric form following Berger, Aftosmis, and Murman 2005 30688f5f89eSMatthew G. Knepley $ The classical flux-limited formulation is psi(r) where 30788f5f89eSMatthew G. Knepley $ 30888f5f89eSMatthew G. Knepley $ r = (u[0] - u[-1]) / (u[1] - u[0]) 30988f5f89eSMatthew G. Knepley $ 31088f5f89eSMatthew G. Knepley $ The second order TVD region is bounded by 31188f5f89eSMatthew G. Knepley $ 31288f5f89eSMatthew G. Knepley $ psi_minmod(r) = min(r,1) and psi_superbee(r) = min(2, 2r, max(1,r)) 31388f5f89eSMatthew G. Knepley $ 31488f5f89eSMatthew G. Knepley $ where all limiters are implicitly clipped to be non-negative. A more convenient slope-limited form is psi(r) = 31588f5f89eSMatthew G. Knepley $ phi(r)(r+1)/2 in which the reconstructed interface values are 31688f5f89eSMatthew G. Knepley $ 31788f5f89eSMatthew G. Knepley $ u(v) = u[0] + phi(r) (grad u)[0] v 31888f5f89eSMatthew G. Knepley $ 31988f5f89eSMatthew G. Knepley $ where v is the vector from centroid to quadrature point. In these variables, the usual limiters become 32088f5f89eSMatthew G. Knepley $ 32188f5f89eSMatthew G. Knepley $ phi_minmod(r) = 2 min(1/(1+r),r/(1+r)) phi_superbee(r) = 2 min(2/(1+r), 2r/(1+r), max(1,r)/(1+r)) 32288f5f89eSMatthew G. Knepley $ 32388f5f89eSMatthew G. Knepley $ For a nicer symmetric formulation, rewrite in terms of 32488f5f89eSMatthew G. Knepley $ 32588f5f89eSMatthew G. Knepley $ f = (u[0] - u[-1]) / (u[1] - u[-1]) 32688f5f89eSMatthew G. Knepley $ 32788f5f89eSMatthew G. Knepley $ where r(f) = f/(1-f). Not that r(1-f) = (1-f)/f = 1/r(f) so the symmetry condition 32888f5f89eSMatthew G. Knepley $ 32988f5f89eSMatthew G. Knepley $ phi(r) = phi(1/r) 33088f5f89eSMatthew G. Knepley $ 33188f5f89eSMatthew G. Knepley $ becomes 33288f5f89eSMatthew G. Knepley $ 33388f5f89eSMatthew G. Knepley $ w(f) = w(1-f). 33488f5f89eSMatthew G. Knepley $ 33588f5f89eSMatthew G. Knepley $ The limiters below implement this final form w(f). The reference methods are 33688f5f89eSMatthew G. Knepley $ 33788f5f89eSMatthew G. Knepley $ w_minmod(f) = 2 min(f,(1-f)) w_superbee(r) = 4 min((1-f), f) 33888f5f89eSMatthew G. Knepley 33988f5f89eSMatthew G. Knepley Level: beginner 34088f5f89eSMatthew G. Knepley 34188f5f89eSMatthew G. Knepley .seealso: PetscLimiterSetType(), PetscLimiterCreate() 34288f5f89eSMatthew G. Knepley @*/ 343c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterLimit(PetscLimiter lim, PetscReal flim, PetscReal *phi) 344c0ce001eSMatthew G. Knepley { 345c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 346c0ce001eSMatthew G. Knepley 347c0ce001eSMatthew G. Knepley PetscFunctionBegin; 348c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 349c0ce001eSMatthew G. Knepley PetscValidPointer(phi, 3); 350c0ce001eSMatthew G. Knepley ierr = (*lim->ops->limit)(lim, flim, phi);CHKERRQ(ierr); 351c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 352c0ce001eSMatthew G. Knepley } 353c0ce001eSMatthew G. Knepley 35488f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterDestroy_Sin(PetscLimiter lim) 355c0ce001eSMatthew G. Knepley { 356c0ce001eSMatthew G. Knepley PetscLimiter_Sin *l = (PetscLimiter_Sin *) lim->data; 357c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 358c0ce001eSMatthew G. Knepley 359c0ce001eSMatthew G. Knepley PetscFunctionBegin; 360c0ce001eSMatthew G. Knepley ierr = PetscFree(l);CHKERRQ(ierr); 361c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 362c0ce001eSMatthew G. Knepley } 363c0ce001eSMatthew G. Knepley 36488f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_Sin_Ascii(PetscLimiter lim, PetscViewer viewer) 365c0ce001eSMatthew G. Knepley { 366c0ce001eSMatthew G. Knepley PetscViewerFormat format; 367c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 368c0ce001eSMatthew G. Knepley 369c0ce001eSMatthew G. Knepley PetscFunctionBegin; 370c0ce001eSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 371c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Sin Slope Limiter:\n");CHKERRQ(ierr); 372c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 373c0ce001eSMatthew G. Knepley } 374c0ce001eSMatthew G. Knepley 37588f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_Sin(PetscLimiter lim, PetscViewer viewer) 376c0ce001eSMatthew G. Knepley { 377c0ce001eSMatthew G. Knepley PetscBool iascii; 378c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 379c0ce001eSMatthew G. Knepley 380c0ce001eSMatthew G. Knepley PetscFunctionBegin; 381c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 382c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 383c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 384c0ce001eSMatthew G. Knepley if (iascii) {ierr = PetscLimiterView_Sin_Ascii(lim, viewer);CHKERRQ(ierr);} 385c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 386c0ce001eSMatthew G. Knepley } 387c0ce001eSMatthew G. Knepley 38888f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterLimit_Sin(PetscLimiter lim, PetscReal f, PetscReal *phi) 389c0ce001eSMatthew G. Knepley { 390c0ce001eSMatthew G. Knepley PetscFunctionBegin; 391c0ce001eSMatthew G. Knepley *phi = PetscSinReal(PETSC_PI*PetscMax(0, PetscMin(f, 1))); 392c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 393c0ce001eSMatthew G. Knepley } 394c0ce001eSMatthew G. Knepley 39588f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterInitialize_Sin(PetscLimiter lim) 396c0ce001eSMatthew G. Knepley { 397c0ce001eSMatthew G. Knepley PetscFunctionBegin; 398c0ce001eSMatthew G. Knepley lim->ops->view = PetscLimiterView_Sin; 399c0ce001eSMatthew G. Knepley lim->ops->destroy = PetscLimiterDestroy_Sin; 400c0ce001eSMatthew G. Knepley lim->ops->limit = PetscLimiterLimit_Sin; 401c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 402c0ce001eSMatthew G. Knepley } 403c0ce001eSMatthew G. Knepley 404c0ce001eSMatthew G. Knepley /*MC 405c0ce001eSMatthew G. Knepley PETSCLIMITERSIN = "sin" - A PetscLimiter object 406c0ce001eSMatthew G. Knepley 407c0ce001eSMatthew G. Knepley Level: intermediate 408c0ce001eSMatthew G. Knepley 409c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 410c0ce001eSMatthew G. Knepley M*/ 411c0ce001eSMatthew G. Knepley 412c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Sin(PetscLimiter lim) 413c0ce001eSMatthew G. Knepley { 414c0ce001eSMatthew G. Knepley PetscLimiter_Sin *l; 415c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 416c0ce001eSMatthew G. Knepley 417c0ce001eSMatthew G. Knepley PetscFunctionBegin; 418c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 419c0ce001eSMatthew G. Knepley ierr = PetscNewLog(lim, &l);CHKERRQ(ierr); 420c0ce001eSMatthew G. Knepley lim->data = l; 421c0ce001eSMatthew G. Knepley 422c0ce001eSMatthew G. Knepley ierr = PetscLimiterInitialize_Sin(lim);CHKERRQ(ierr); 423c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 424c0ce001eSMatthew G. Knepley } 425c0ce001eSMatthew G. Knepley 42688f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterDestroy_Zero(PetscLimiter lim) 427c0ce001eSMatthew G. Knepley { 428c0ce001eSMatthew G. Knepley PetscLimiter_Zero *l = (PetscLimiter_Zero *) lim->data; 429c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 430c0ce001eSMatthew G. Knepley 431c0ce001eSMatthew G. Knepley PetscFunctionBegin; 432c0ce001eSMatthew G. Knepley ierr = PetscFree(l);CHKERRQ(ierr); 433c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 434c0ce001eSMatthew G. Knepley } 435c0ce001eSMatthew G. Knepley 43688f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_Zero_Ascii(PetscLimiter lim, PetscViewer viewer) 437c0ce001eSMatthew G. Knepley { 438c0ce001eSMatthew G. Knepley PetscViewerFormat format; 439c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 440c0ce001eSMatthew G. Knepley 441c0ce001eSMatthew G. Knepley PetscFunctionBegin; 442c0ce001eSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 443c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Zero Slope Limiter:\n");CHKERRQ(ierr); 444c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 445c0ce001eSMatthew G. Knepley } 446c0ce001eSMatthew G. Knepley 44788f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_Zero(PetscLimiter lim, PetscViewer viewer) 448c0ce001eSMatthew G. Knepley { 449c0ce001eSMatthew G. Knepley PetscBool iascii; 450c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 451c0ce001eSMatthew G. Knepley 452c0ce001eSMatthew G. Knepley PetscFunctionBegin; 453c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 454c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 455c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 456c0ce001eSMatthew G. Knepley if (iascii) {ierr = PetscLimiterView_Zero_Ascii(lim, viewer);CHKERRQ(ierr);} 457c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 458c0ce001eSMatthew G. Knepley } 459c0ce001eSMatthew G. Knepley 46088f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterLimit_Zero(PetscLimiter lim, PetscReal f, PetscReal *phi) 461c0ce001eSMatthew G. Knepley { 462c0ce001eSMatthew G. Knepley PetscFunctionBegin; 463c0ce001eSMatthew G. Knepley *phi = 0.0; 464c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 465c0ce001eSMatthew G. Knepley } 466c0ce001eSMatthew G. Knepley 46788f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterInitialize_Zero(PetscLimiter lim) 468c0ce001eSMatthew G. Knepley { 469c0ce001eSMatthew G. Knepley PetscFunctionBegin; 470c0ce001eSMatthew G. Knepley lim->ops->view = PetscLimiterView_Zero; 471c0ce001eSMatthew G. Knepley lim->ops->destroy = PetscLimiterDestroy_Zero; 472c0ce001eSMatthew G. Knepley lim->ops->limit = PetscLimiterLimit_Zero; 473c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 474c0ce001eSMatthew G. Knepley } 475c0ce001eSMatthew G. Knepley 476c0ce001eSMatthew G. Knepley /*MC 477c0ce001eSMatthew G. Knepley PETSCLIMITERZERO = "zero" - A PetscLimiter object 478c0ce001eSMatthew G. Knepley 479c0ce001eSMatthew G. Knepley Level: intermediate 480c0ce001eSMatthew G. Knepley 481c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 482c0ce001eSMatthew G. Knepley M*/ 483c0ce001eSMatthew G. Knepley 484c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Zero(PetscLimiter lim) 485c0ce001eSMatthew G. Knepley { 486c0ce001eSMatthew G. Knepley PetscLimiter_Zero *l; 487c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 488c0ce001eSMatthew G. Knepley 489c0ce001eSMatthew G. Knepley PetscFunctionBegin; 490c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 491c0ce001eSMatthew G. Knepley ierr = PetscNewLog(lim, &l);CHKERRQ(ierr); 492c0ce001eSMatthew G. Knepley lim->data = l; 493c0ce001eSMatthew G. Knepley 494c0ce001eSMatthew G. Knepley ierr = PetscLimiterInitialize_Zero(lim);CHKERRQ(ierr); 495c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 496c0ce001eSMatthew G. Knepley } 497c0ce001eSMatthew G. Knepley 49888f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterDestroy_None(PetscLimiter lim) 499c0ce001eSMatthew G. Knepley { 500c0ce001eSMatthew G. Knepley PetscLimiter_None *l = (PetscLimiter_None *) lim->data; 501c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 502c0ce001eSMatthew G. Knepley 503c0ce001eSMatthew G. Knepley PetscFunctionBegin; 504c0ce001eSMatthew G. Knepley ierr = PetscFree(l);CHKERRQ(ierr); 505c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 506c0ce001eSMatthew G. Knepley } 507c0ce001eSMatthew G. Knepley 50888f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_None_Ascii(PetscLimiter lim, PetscViewer viewer) 509c0ce001eSMatthew G. Knepley { 510c0ce001eSMatthew G. Knepley PetscViewerFormat format; 511c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 512c0ce001eSMatthew G. Knepley 513c0ce001eSMatthew G. Knepley PetscFunctionBegin; 514c0ce001eSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 515c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "None Slope Limiter:\n");CHKERRQ(ierr); 516c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 517c0ce001eSMatthew G. Knepley } 518c0ce001eSMatthew G. Knepley 51988f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_None(PetscLimiter lim, PetscViewer viewer) 520c0ce001eSMatthew G. Knepley { 521c0ce001eSMatthew G. Knepley PetscBool iascii; 522c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 523c0ce001eSMatthew G. Knepley 524c0ce001eSMatthew G. Knepley PetscFunctionBegin; 525c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 526c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 527c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 528c0ce001eSMatthew G. Knepley if (iascii) {ierr = PetscLimiterView_None_Ascii(lim, viewer);CHKERRQ(ierr);} 529c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 530c0ce001eSMatthew G. Knepley } 531c0ce001eSMatthew G. Knepley 53288f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterLimit_None(PetscLimiter lim, PetscReal f, PetscReal *phi) 533c0ce001eSMatthew G. Knepley { 534c0ce001eSMatthew G. Knepley PetscFunctionBegin; 535c0ce001eSMatthew G. Knepley *phi = 1.0; 536c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 537c0ce001eSMatthew G. Knepley } 538c0ce001eSMatthew G. Knepley 53988f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterInitialize_None(PetscLimiter lim) 540c0ce001eSMatthew G. Knepley { 541c0ce001eSMatthew G. Knepley PetscFunctionBegin; 542c0ce001eSMatthew G. Knepley lim->ops->view = PetscLimiterView_None; 543c0ce001eSMatthew G. Knepley lim->ops->destroy = PetscLimiterDestroy_None; 544c0ce001eSMatthew G. Knepley lim->ops->limit = PetscLimiterLimit_None; 545c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 546c0ce001eSMatthew G. Knepley } 547c0ce001eSMatthew G. Knepley 548c0ce001eSMatthew G. Knepley /*MC 549c0ce001eSMatthew G. Knepley PETSCLIMITERNONE = "none" - A PetscLimiter object 550c0ce001eSMatthew G. Knepley 551c0ce001eSMatthew G. Knepley Level: intermediate 552c0ce001eSMatthew G. Knepley 553c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 554c0ce001eSMatthew G. Knepley M*/ 555c0ce001eSMatthew G. Knepley 556c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_None(PetscLimiter lim) 557c0ce001eSMatthew G. Knepley { 558c0ce001eSMatthew G. Knepley PetscLimiter_None *l; 559c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 560c0ce001eSMatthew G. Knepley 561c0ce001eSMatthew G. Knepley PetscFunctionBegin; 562c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 563c0ce001eSMatthew G. Knepley ierr = PetscNewLog(lim, &l);CHKERRQ(ierr); 564c0ce001eSMatthew G. Knepley lim->data = l; 565c0ce001eSMatthew G. Knepley 566c0ce001eSMatthew G. Knepley ierr = PetscLimiterInitialize_None(lim);CHKERRQ(ierr); 567c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 568c0ce001eSMatthew G. Knepley } 569c0ce001eSMatthew G. Knepley 57088f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterDestroy_Minmod(PetscLimiter lim) 571c0ce001eSMatthew G. Knepley { 572c0ce001eSMatthew G. Knepley PetscLimiter_Minmod *l = (PetscLimiter_Minmod *) lim->data; 573c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 574c0ce001eSMatthew G. Knepley 575c0ce001eSMatthew G. Knepley PetscFunctionBegin; 576c0ce001eSMatthew G. Knepley ierr = PetscFree(l);CHKERRQ(ierr); 577c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 578c0ce001eSMatthew G. Knepley } 579c0ce001eSMatthew G. Knepley 58088f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_Minmod_Ascii(PetscLimiter lim, PetscViewer viewer) 581c0ce001eSMatthew G. Knepley { 582c0ce001eSMatthew G. Knepley PetscViewerFormat format; 583c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 584c0ce001eSMatthew G. Knepley 585c0ce001eSMatthew G. Knepley PetscFunctionBegin; 586c0ce001eSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 587c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Minmod Slope Limiter:\n");CHKERRQ(ierr); 588c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 589c0ce001eSMatthew G. Knepley } 590c0ce001eSMatthew G. Knepley 59188f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_Minmod(PetscLimiter lim, PetscViewer viewer) 592c0ce001eSMatthew G. Knepley { 593c0ce001eSMatthew G. Knepley PetscBool iascii; 594c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 595c0ce001eSMatthew G. Knepley 596c0ce001eSMatthew G. Knepley PetscFunctionBegin; 597c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 598c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 599c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 600c0ce001eSMatthew G. Knepley if (iascii) {ierr = PetscLimiterView_Minmod_Ascii(lim, viewer);CHKERRQ(ierr);} 601c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 602c0ce001eSMatthew G. Knepley } 603c0ce001eSMatthew G. Knepley 60488f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterLimit_Minmod(PetscLimiter lim, PetscReal f, PetscReal *phi) 605c0ce001eSMatthew G. Knepley { 606c0ce001eSMatthew G. Knepley PetscFunctionBegin; 607c0ce001eSMatthew G. Knepley *phi = 2*PetscMax(0, PetscMin(f, 1-f)); 608c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 609c0ce001eSMatthew G. Knepley } 610c0ce001eSMatthew G. Knepley 61188f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterInitialize_Minmod(PetscLimiter lim) 612c0ce001eSMatthew G. Knepley { 613c0ce001eSMatthew G. Knepley PetscFunctionBegin; 614c0ce001eSMatthew G. Knepley lim->ops->view = PetscLimiterView_Minmod; 615c0ce001eSMatthew G. Knepley lim->ops->destroy = PetscLimiterDestroy_Minmod; 616c0ce001eSMatthew G. Knepley lim->ops->limit = PetscLimiterLimit_Minmod; 617c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 618c0ce001eSMatthew G. Knepley } 619c0ce001eSMatthew G. Knepley 620c0ce001eSMatthew G. Knepley /*MC 621c0ce001eSMatthew G. Knepley PETSCLIMITERMINMOD = "minmod" - A PetscLimiter object 622c0ce001eSMatthew G. Knepley 623c0ce001eSMatthew G. Knepley Level: intermediate 624c0ce001eSMatthew G. Knepley 625c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 626c0ce001eSMatthew G. Knepley M*/ 627c0ce001eSMatthew G. Knepley 628c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Minmod(PetscLimiter lim) 629c0ce001eSMatthew G. Knepley { 630c0ce001eSMatthew G. Knepley PetscLimiter_Minmod *l; 631c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 632c0ce001eSMatthew G. Knepley 633c0ce001eSMatthew G. Knepley PetscFunctionBegin; 634c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 635c0ce001eSMatthew G. Knepley ierr = PetscNewLog(lim, &l);CHKERRQ(ierr); 636c0ce001eSMatthew G. Knepley lim->data = l; 637c0ce001eSMatthew G. Knepley 638c0ce001eSMatthew G. Knepley ierr = PetscLimiterInitialize_Minmod(lim);CHKERRQ(ierr); 639c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 640c0ce001eSMatthew G. Knepley } 641c0ce001eSMatthew G. Knepley 64288f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterDestroy_VanLeer(PetscLimiter lim) 643c0ce001eSMatthew G. Knepley { 644c0ce001eSMatthew G. Knepley PetscLimiter_VanLeer *l = (PetscLimiter_VanLeer *) lim->data; 645c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 646c0ce001eSMatthew G. Knepley 647c0ce001eSMatthew G. Knepley PetscFunctionBegin; 648c0ce001eSMatthew G. Knepley ierr = PetscFree(l);CHKERRQ(ierr); 649c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 650c0ce001eSMatthew G. Knepley } 651c0ce001eSMatthew G. Knepley 65288f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_VanLeer_Ascii(PetscLimiter lim, PetscViewer viewer) 653c0ce001eSMatthew G. Knepley { 654c0ce001eSMatthew G. Knepley PetscViewerFormat format; 655c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 656c0ce001eSMatthew G. Knepley 657c0ce001eSMatthew G. Knepley PetscFunctionBegin; 658c0ce001eSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 659c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Van Leer Slope Limiter:\n");CHKERRQ(ierr); 660c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 661c0ce001eSMatthew G. Knepley } 662c0ce001eSMatthew G. Knepley 66388f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_VanLeer(PetscLimiter lim, PetscViewer viewer) 664c0ce001eSMatthew G. Knepley { 665c0ce001eSMatthew G. Knepley PetscBool iascii; 666c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 667c0ce001eSMatthew G. Knepley 668c0ce001eSMatthew G. Knepley PetscFunctionBegin; 669c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 670c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 671c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 672c0ce001eSMatthew G. Knepley if (iascii) {ierr = PetscLimiterView_VanLeer_Ascii(lim, viewer);CHKERRQ(ierr);} 673c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 674c0ce001eSMatthew G. Knepley } 675c0ce001eSMatthew G. Knepley 67688f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterLimit_VanLeer(PetscLimiter lim, PetscReal f, PetscReal *phi) 677c0ce001eSMatthew G. Knepley { 678c0ce001eSMatthew G. Knepley PetscFunctionBegin; 679c0ce001eSMatthew G. Knepley *phi = PetscMax(0, 4*f*(1-f)); 680c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 681c0ce001eSMatthew G. Knepley } 682c0ce001eSMatthew G. Knepley 68388f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterInitialize_VanLeer(PetscLimiter lim) 684c0ce001eSMatthew G. Knepley { 685c0ce001eSMatthew G. Knepley PetscFunctionBegin; 686c0ce001eSMatthew G. Knepley lim->ops->view = PetscLimiterView_VanLeer; 687c0ce001eSMatthew G. Knepley lim->ops->destroy = PetscLimiterDestroy_VanLeer; 688c0ce001eSMatthew G. Knepley lim->ops->limit = PetscLimiterLimit_VanLeer; 689c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 690c0ce001eSMatthew G. Knepley } 691c0ce001eSMatthew G. Knepley 692c0ce001eSMatthew G. Knepley /*MC 693c0ce001eSMatthew G. Knepley PETSCLIMITERVANLEER = "vanleer" - A PetscLimiter object 694c0ce001eSMatthew G. Knepley 695c0ce001eSMatthew G. Knepley Level: intermediate 696c0ce001eSMatthew G. Knepley 697c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 698c0ce001eSMatthew G. Knepley M*/ 699c0ce001eSMatthew G. Knepley 700c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_VanLeer(PetscLimiter lim) 701c0ce001eSMatthew G. Knepley { 702c0ce001eSMatthew G. Knepley PetscLimiter_VanLeer *l; 703c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 704c0ce001eSMatthew G. Knepley 705c0ce001eSMatthew G. Knepley PetscFunctionBegin; 706c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 707c0ce001eSMatthew G. Knepley ierr = PetscNewLog(lim, &l);CHKERRQ(ierr); 708c0ce001eSMatthew G. Knepley lim->data = l; 709c0ce001eSMatthew G. Knepley 710c0ce001eSMatthew G. Knepley ierr = PetscLimiterInitialize_VanLeer(lim);CHKERRQ(ierr); 711c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 712c0ce001eSMatthew G. Knepley } 713c0ce001eSMatthew G. Knepley 71488f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterDestroy_VanAlbada(PetscLimiter lim) 715c0ce001eSMatthew G. Knepley { 716c0ce001eSMatthew G. Knepley PetscLimiter_VanAlbada *l = (PetscLimiter_VanAlbada *) lim->data; 717c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 718c0ce001eSMatthew G. Knepley 719c0ce001eSMatthew G. Knepley PetscFunctionBegin; 720c0ce001eSMatthew G. Knepley ierr = PetscFree(l);CHKERRQ(ierr); 721c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 722c0ce001eSMatthew G. Knepley } 723c0ce001eSMatthew G. Knepley 72488f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_VanAlbada_Ascii(PetscLimiter lim, PetscViewer viewer) 725c0ce001eSMatthew G. Knepley { 726c0ce001eSMatthew G. Knepley PetscViewerFormat format; 727c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 728c0ce001eSMatthew G. Knepley 729c0ce001eSMatthew G. Knepley PetscFunctionBegin; 730c0ce001eSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 731c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Van Albada Slope Limiter:\n");CHKERRQ(ierr); 732c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 733c0ce001eSMatthew G. Knepley } 734c0ce001eSMatthew G. Knepley 73588f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_VanAlbada(PetscLimiter lim, PetscViewer viewer) 736c0ce001eSMatthew G. Knepley { 737c0ce001eSMatthew G. Knepley PetscBool iascii; 738c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 739c0ce001eSMatthew G. Knepley 740c0ce001eSMatthew G. Knepley PetscFunctionBegin; 741c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 742c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 743c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 744c0ce001eSMatthew G. Knepley if (iascii) {ierr = PetscLimiterView_VanAlbada_Ascii(lim, viewer);CHKERRQ(ierr);} 745c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 746c0ce001eSMatthew G. Knepley } 747c0ce001eSMatthew G. Knepley 74888f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterLimit_VanAlbada(PetscLimiter lim, PetscReal f, PetscReal *phi) 749c0ce001eSMatthew G. Knepley { 750c0ce001eSMatthew G. Knepley PetscFunctionBegin; 751c0ce001eSMatthew G. Knepley *phi = PetscMax(0, 2*f*(1-f) / (PetscSqr(f) + PetscSqr(1-f))); 752c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 753c0ce001eSMatthew G. Knepley } 754c0ce001eSMatthew G. Knepley 75588f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterInitialize_VanAlbada(PetscLimiter lim) 756c0ce001eSMatthew G. Knepley { 757c0ce001eSMatthew G. Knepley PetscFunctionBegin; 758c0ce001eSMatthew G. Knepley lim->ops->view = PetscLimiterView_VanAlbada; 759c0ce001eSMatthew G. Knepley lim->ops->destroy = PetscLimiterDestroy_VanAlbada; 760c0ce001eSMatthew G. Knepley lim->ops->limit = PetscLimiterLimit_VanAlbada; 761c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 762c0ce001eSMatthew G. Knepley } 763c0ce001eSMatthew G. Knepley 764c0ce001eSMatthew G. Knepley /*MC 765c0ce001eSMatthew G. Knepley PETSCLIMITERVANALBADA = "vanalbada" - A PetscLimiter object 766c0ce001eSMatthew G. Knepley 767c0ce001eSMatthew G. Knepley Level: intermediate 768c0ce001eSMatthew G. Knepley 769c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 770c0ce001eSMatthew G. Knepley M*/ 771c0ce001eSMatthew G. Knepley 772c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_VanAlbada(PetscLimiter lim) 773c0ce001eSMatthew G. Knepley { 774c0ce001eSMatthew G. Knepley PetscLimiter_VanAlbada *l; 775c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 776c0ce001eSMatthew G. Knepley 777c0ce001eSMatthew G. Knepley PetscFunctionBegin; 778c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 779c0ce001eSMatthew G. Knepley ierr = PetscNewLog(lim, &l);CHKERRQ(ierr); 780c0ce001eSMatthew G. Knepley lim->data = l; 781c0ce001eSMatthew G. Knepley 782c0ce001eSMatthew G. Knepley ierr = PetscLimiterInitialize_VanAlbada(lim);CHKERRQ(ierr); 783c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 784c0ce001eSMatthew G. Knepley } 785c0ce001eSMatthew G. Knepley 78688f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterDestroy_Superbee(PetscLimiter lim) 787c0ce001eSMatthew G. Knepley { 788c0ce001eSMatthew G. Knepley PetscLimiter_Superbee *l = (PetscLimiter_Superbee *) lim->data; 789c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 790c0ce001eSMatthew G. Knepley 791c0ce001eSMatthew G. Knepley PetscFunctionBegin; 792c0ce001eSMatthew G. Knepley ierr = PetscFree(l);CHKERRQ(ierr); 793c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 794c0ce001eSMatthew G. Knepley } 795c0ce001eSMatthew G. Knepley 79688f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_Superbee_Ascii(PetscLimiter lim, PetscViewer viewer) 797c0ce001eSMatthew G. Knepley { 798c0ce001eSMatthew G. Knepley PetscViewerFormat format; 799c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 800c0ce001eSMatthew G. Knepley 801c0ce001eSMatthew G. Knepley PetscFunctionBegin; 802c0ce001eSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 803c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Superbee Slope Limiter:\n");CHKERRQ(ierr); 804c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 805c0ce001eSMatthew G. Knepley } 806c0ce001eSMatthew G. Knepley 80788f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_Superbee(PetscLimiter lim, PetscViewer viewer) 808c0ce001eSMatthew G. Knepley { 809c0ce001eSMatthew G. Knepley PetscBool iascii; 810c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 811c0ce001eSMatthew G. Knepley 812c0ce001eSMatthew G. Knepley PetscFunctionBegin; 813c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 814c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 815c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 816c0ce001eSMatthew G. Knepley if (iascii) {ierr = PetscLimiterView_Superbee_Ascii(lim, viewer);CHKERRQ(ierr);} 817c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 818c0ce001eSMatthew G. Knepley } 819c0ce001eSMatthew G. Knepley 82088f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterLimit_Superbee(PetscLimiter lim, PetscReal f, PetscReal *phi) 821c0ce001eSMatthew G. Knepley { 822c0ce001eSMatthew G. Knepley PetscFunctionBegin; 823c0ce001eSMatthew G. Knepley *phi = 4*PetscMax(0, PetscMin(f, 1-f)); 824c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 825c0ce001eSMatthew G. Knepley } 826c0ce001eSMatthew G. Knepley 82788f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterInitialize_Superbee(PetscLimiter lim) 828c0ce001eSMatthew G. Knepley { 829c0ce001eSMatthew G. Knepley PetscFunctionBegin; 830c0ce001eSMatthew G. Knepley lim->ops->view = PetscLimiterView_Superbee; 831c0ce001eSMatthew G. Knepley lim->ops->destroy = PetscLimiterDestroy_Superbee; 832c0ce001eSMatthew G. Knepley lim->ops->limit = PetscLimiterLimit_Superbee; 833c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 834c0ce001eSMatthew G. Knepley } 835c0ce001eSMatthew G. Knepley 836c0ce001eSMatthew G. Knepley /*MC 837c0ce001eSMatthew G. Knepley PETSCLIMITERSUPERBEE = "superbee" - A PetscLimiter object 838c0ce001eSMatthew G. Knepley 839c0ce001eSMatthew G. Knepley Level: intermediate 840c0ce001eSMatthew G. Knepley 841c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 842c0ce001eSMatthew G. Knepley M*/ 843c0ce001eSMatthew G. Knepley 844c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Superbee(PetscLimiter lim) 845c0ce001eSMatthew G. Knepley { 846c0ce001eSMatthew G. Knepley PetscLimiter_Superbee *l; 847c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 848c0ce001eSMatthew G. Knepley 849c0ce001eSMatthew G. Knepley PetscFunctionBegin; 850c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 851c0ce001eSMatthew G. Knepley ierr = PetscNewLog(lim, &l);CHKERRQ(ierr); 852c0ce001eSMatthew G. Knepley lim->data = l; 853c0ce001eSMatthew G. Knepley 854c0ce001eSMatthew G. Knepley ierr = PetscLimiterInitialize_Superbee(lim);CHKERRQ(ierr); 855c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 856c0ce001eSMatthew G. Knepley } 857c0ce001eSMatthew G. Knepley 85888f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterDestroy_MC(PetscLimiter lim) 859c0ce001eSMatthew G. Knepley { 860c0ce001eSMatthew G. Knepley PetscLimiter_MC *l = (PetscLimiter_MC *) lim->data; 861c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 862c0ce001eSMatthew G. Knepley 863c0ce001eSMatthew G. Knepley PetscFunctionBegin; 864c0ce001eSMatthew G. Knepley ierr = PetscFree(l);CHKERRQ(ierr); 865c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 866c0ce001eSMatthew G. Knepley } 867c0ce001eSMatthew G. Knepley 86888f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_MC_Ascii(PetscLimiter lim, PetscViewer viewer) 869c0ce001eSMatthew G. Knepley { 870c0ce001eSMatthew G. Knepley PetscViewerFormat format; 871c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 872c0ce001eSMatthew G. Knepley 873c0ce001eSMatthew G. Knepley PetscFunctionBegin; 874c0ce001eSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 875c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "MC Slope Limiter:\n");CHKERRQ(ierr); 876c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 877c0ce001eSMatthew G. Knepley } 878c0ce001eSMatthew G. Knepley 87988f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterView_MC(PetscLimiter lim, PetscViewer viewer) 880c0ce001eSMatthew G. Knepley { 881c0ce001eSMatthew G. Knepley PetscBool iascii; 882c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 883c0ce001eSMatthew G. Knepley 884c0ce001eSMatthew G. Knepley PetscFunctionBegin; 885c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 886c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 887c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 888c0ce001eSMatthew G. Knepley if (iascii) {ierr = PetscLimiterView_MC_Ascii(lim, viewer);CHKERRQ(ierr);} 889c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 890c0ce001eSMatthew G. Knepley } 891c0ce001eSMatthew G. Knepley 892c0ce001eSMatthew G. Knepley /* aka Barth-Jespersen */ 89388f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterLimit_MC(PetscLimiter lim, PetscReal f, PetscReal *phi) 894c0ce001eSMatthew G. Knepley { 895c0ce001eSMatthew G. Knepley PetscFunctionBegin; 896c0ce001eSMatthew G. Knepley *phi = PetscMin(1, 4*PetscMax(0, PetscMin(f, 1-f))); 897c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 898c0ce001eSMatthew G. Knepley } 899c0ce001eSMatthew G. Knepley 90088f5f89eSMatthew G. Knepley static PetscErrorCode PetscLimiterInitialize_MC(PetscLimiter lim) 901c0ce001eSMatthew G. Knepley { 902c0ce001eSMatthew G. Knepley PetscFunctionBegin; 903c0ce001eSMatthew G. Knepley lim->ops->view = PetscLimiterView_MC; 904c0ce001eSMatthew G. Knepley lim->ops->destroy = PetscLimiterDestroy_MC; 905c0ce001eSMatthew G. Knepley lim->ops->limit = PetscLimiterLimit_MC; 906c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 907c0ce001eSMatthew G. Knepley } 908c0ce001eSMatthew G. Knepley 909c0ce001eSMatthew G. Knepley /*MC 910c0ce001eSMatthew G. Knepley PETSCLIMITERMC = "mc" - A PetscLimiter object 911c0ce001eSMatthew G. Knepley 912c0ce001eSMatthew G. Knepley Level: intermediate 913c0ce001eSMatthew G. Knepley 914c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 915c0ce001eSMatthew G. Knepley M*/ 916c0ce001eSMatthew G. Knepley 917c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_MC(PetscLimiter lim) 918c0ce001eSMatthew G. Knepley { 919c0ce001eSMatthew G. Knepley PetscLimiter_MC *l; 920c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 921c0ce001eSMatthew G. Knepley 922c0ce001eSMatthew G. Knepley PetscFunctionBegin; 923c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 924c0ce001eSMatthew G. Knepley ierr = PetscNewLog(lim, &l);CHKERRQ(ierr); 925c0ce001eSMatthew G. Knepley lim->data = l; 926c0ce001eSMatthew G. Knepley 927c0ce001eSMatthew G. Knepley ierr = PetscLimiterInitialize_MC(lim);CHKERRQ(ierr); 928c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 929c0ce001eSMatthew G. Knepley } 930c0ce001eSMatthew G. Knepley 931c0ce001eSMatthew G. Knepley PetscClassId PETSCFV_CLASSID = 0; 932c0ce001eSMatthew G. Knepley 933c0ce001eSMatthew G. Knepley PetscFunctionList PetscFVList = NULL; 934c0ce001eSMatthew G. Knepley PetscBool PetscFVRegisterAllCalled = PETSC_FALSE; 935c0ce001eSMatthew G. Knepley 936c0ce001eSMatthew G. Knepley /*@C 937c0ce001eSMatthew G. Knepley PetscFVRegister - Adds a new PetscFV implementation 938c0ce001eSMatthew G. Knepley 939c0ce001eSMatthew G. Knepley Not Collective 940c0ce001eSMatthew G. Knepley 941c0ce001eSMatthew G. Knepley Input Parameters: 942c0ce001eSMatthew G. Knepley + name - The name of a new user-defined creation routine 943c0ce001eSMatthew G. Knepley - create_func - The creation routine itself 944c0ce001eSMatthew G. Knepley 945c0ce001eSMatthew G. Knepley Notes: 946c0ce001eSMatthew G. Knepley PetscFVRegister() may be called multiple times to add several user-defined PetscFVs 947c0ce001eSMatthew G. Knepley 948c0ce001eSMatthew G. Knepley Sample usage: 949c0ce001eSMatthew G. Knepley .vb 950c0ce001eSMatthew G. Knepley PetscFVRegister("my_fv", MyPetscFVCreate); 951c0ce001eSMatthew G. Knepley .ve 952c0ce001eSMatthew G. Knepley 953c0ce001eSMatthew G. Knepley Then, your PetscFV type can be chosen with the procedural interface via 954c0ce001eSMatthew G. Knepley .vb 955c0ce001eSMatthew G. Knepley PetscFVCreate(MPI_Comm, PetscFV *); 956c0ce001eSMatthew G. Knepley PetscFVSetType(PetscFV, "my_fv"); 957c0ce001eSMatthew G. Knepley .ve 958c0ce001eSMatthew G. Knepley or at runtime via the option 959c0ce001eSMatthew G. Knepley .vb 960c0ce001eSMatthew G. Knepley -petscfv_type my_fv 961c0ce001eSMatthew G. Knepley .ve 962c0ce001eSMatthew G. Knepley 963c0ce001eSMatthew G. Knepley Level: advanced 964c0ce001eSMatthew G. Knepley 965c0ce001eSMatthew G. Knepley .seealso: PetscFVRegisterAll(), PetscFVRegisterDestroy() 966c0ce001eSMatthew G. Knepley 967c0ce001eSMatthew G. Knepley @*/ 968c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVRegister(const char sname[], PetscErrorCode (*function)(PetscFV)) 969c0ce001eSMatthew G. Knepley { 970c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 971c0ce001eSMatthew G. Knepley 972c0ce001eSMatthew G. Knepley PetscFunctionBegin; 973c0ce001eSMatthew G. Knepley ierr = PetscFunctionListAdd(&PetscFVList, sname, function);CHKERRQ(ierr); 974c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 975c0ce001eSMatthew G. Knepley } 976c0ce001eSMatthew G. Knepley 977c0ce001eSMatthew G. Knepley /*@C 978c0ce001eSMatthew G. Knepley PetscFVSetType - Builds a particular PetscFV 979c0ce001eSMatthew G. Knepley 980c0ce001eSMatthew G. Knepley Collective on fvm 981c0ce001eSMatthew G. Knepley 982c0ce001eSMatthew G. Knepley Input Parameters: 983c0ce001eSMatthew G. Knepley + fvm - The PetscFV object 984c0ce001eSMatthew G. Knepley - name - The kind of FVM space 985c0ce001eSMatthew G. Knepley 986c0ce001eSMatthew G. Knepley Options Database Key: 987c0ce001eSMatthew G. Knepley . -petscfv_type <type> - Sets the PetscFV type; use -help for a list of available types 988c0ce001eSMatthew G. Knepley 989c0ce001eSMatthew G. Knepley Level: intermediate 990c0ce001eSMatthew G. Knepley 991c0ce001eSMatthew G. Knepley .seealso: PetscFVGetType(), PetscFVCreate() 992c0ce001eSMatthew G. Knepley @*/ 993c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetType(PetscFV fvm, PetscFVType name) 994c0ce001eSMatthew G. Knepley { 995c0ce001eSMatthew G. Knepley PetscErrorCode (*r)(PetscFV); 996c0ce001eSMatthew G. Knepley PetscBool match; 997c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 998c0ce001eSMatthew G. Knepley 999c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1000c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1001c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) fvm, name, &match);CHKERRQ(ierr); 1002c0ce001eSMatthew G. Knepley if (match) PetscFunctionReturn(0); 1003c0ce001eSMatthew G. Knepley 1004c0ce001eSMatthew G. Knepley ierr = PetscFVRegisterAll();CHKERRQ(ierr); 1005c0ce001eSMatthew G. Knepley ierr = PetscFunctionListFind(PetscFVList, name, &r);CHKERRQ(ierr); 1006c0ce001eSMatthew G. Knepley if (!r) SETERRQ1(PetscObjectComm((PetscObject) fvm), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFV type: %s", name); 1007c0ce001eSMatthew G. Knepley 1008c0ce001eSMatthew G. Knepley if (fvm->ops->destroy) { 1009c0ce001eSMatthew G. Knepley ierr = (*fvm->ops->destroy)(fvm);CHKERRQ(ierr); 1010c0ce001eSMatthew G. Knepley fvm->ops->destroy = NULL; 1011c0ce001eSMatthew G. Knepley } 1012c0ce001eSMatthew G. Knepley ierr = (*r)(fvm);CHKERRQ(ierr); 1013c0ce001eSMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) fvm, name);CHKERRQ(ierr); 1014c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1015c0ce001eSMatthew G. Knepley } 1016c0ce001eSMatthew G. Knepley 1017c0ce001eSMatthew G. Knepley /*@C 1018c0ce001eSMatthew G. Knepley PetscFVGetType - Gets the PetscFV type name (as a string) from the object. 1019c0ce001eSMatthew G. Knepley 1020c0ce001eSMatthew G. Knepley Not Collective 1021c0ce001eSMatthew G. Knepley 1022c0ce001eSMatthew G. Knepley Input Parameter: 1023c0ce001eSMatthew G. Knepley . fvm - The PetscFV 1024c0ce001eSMatthew G. Knepley 1025c0ce001eSMatthew G. Knepley Output Parameter: 1026c0ce001eSMatthew G. Knepley . name - The PetscFV type name 1027c0ce001eSMatthew G. Knepley 1028c0ce001eSMatthew G. Knepley Level: intermediate 1029c0ce001eSMatthew G. Knepley 1030c0ce001eSMatthew G. Knepley .seealso: PetscFVSetType(), PetscFVCreate() 1031c0ce001eSMatthew G. Knepley @*/ 1032c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetType(PetscFV fvm, PetscFVType *name) 1033c0ce001eSMatthew G. Knepley { 1034c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1035c0ce001eSMatthew G. Knepley 1036c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1037c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1038c0ce001eSMatthew G. Knepley PetscValidPointer(name, 2); 1039c0ce001eSMatthew G. Knepley ierr = PetscFVRegisterAll();CHKERRQ(ierr); 1040c0ce001eSMatthew G. Knepley *name = ((PetscObject) fvm)->type_name; 1041c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1042c0ce001eSMatthew G. Knepley } 1043c0ce001eSMatthew G. Knepley 1044c0ce001eSMatthew G. Knepley /*@C 1045fe2efc57SMark PetscFVViewFromOptions - View from Options 1046fe2efc57SMark 1047fe2efc57SMark Collective on PetscFV 1048fe2efc57SMark 1049fe2efc57SMark Input Parameters: 1050fe2efc57SMark + A - the PetscFV object 1051736c3998SJose E. Roman . obj - Optional object 1052736c3998SJose E. Roman - name - command line option 1053fe2efc57SMark 1054fe2efc57SMark Level: intermediate 1055fe2efc57SMark .seealso: PetscFV, PetscFVView, PetscObjectViewFromOptions(), PetscFVCreate() 1056fe2efc57SMark @*/ 1057fe2efc57SMark PetscErrorCode PetscFVViewFromOptions(PetscFV A,PetscObject obj,const char name[]) 1058fe2efc57SMark { 1059fe2efc57SMark PetscErrorCode ierr; 1060fe2efc57SMark 1061fe2efc57SMark PetscFunctionBegin; 1062fe2efc57SMark PetscValidHeaderSpecific(A,PETSCFV_CLASSID,1); 1063fe2efc57SMark ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr); 1064fe2efc57SMark PetscFunctionReturn(0); 1065fe2efc57SMark } 1066fe2efc57SMark 1067fe2efc57SMark /*@C 1068c0ce001eSMatthew G. Knepley PetscFVView - Views a PetscFV 1069c0ce001eSMatthew G. Knepley 1070c0ce001eSMatthew G. Knepley Collective on fvm 1071c0ce001eSMatthew G. Knepley 1072c0ce001eSMatthew G. Knepley Input Parameter: 1073c0ce001eSMatthew G. Knepley + fvm - the PetscFV object to view 1074c0ce001eSMatthew G. Knepley - v - the viewer 1075c0ce001eSMatthew G. Knepley 107688f5f89eSMatthew G. Knepley Level: beginner 1077c0ce001eSMatthew G. Knepley 1078c0ce001eSMatthew G. Knepley .seealso: PetscFVDestroy() 1079c0ce001eSMatthew G. Knepley @*/ 1080c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVView(PetscFV fvm, PetscViewer v) 1081c0ce001eSMatthew G. Knepley { 1082c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1083c0ce001eSMatthew G. Knepley 1084c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1085c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1086c0ce001eSMatthew G. Knepley if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fvm), &v);CHKERRQ(ierr);} 1087c0ce001eSMatthew G. Knepley if (fvm->ops->view) {ierr = (*fvm->ops->view)(fvm, v);CHKERRQ(ierr);} 1088c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1089c0ce001eSMatthew G. Knepley } 1090c0ce001eSMatthew G. Knepley 1091c0ce001eSMatthew G. Knepley /*@ 1092c0ce001eSMatthew G. Knepley PetscFVSetFromOptions - sets parameters in a PetscFV from the options database 1093c0ce001eSMatthew G. Knepley 1094c0ce001eSMatthew G. Knepley Collective on fvm 1095c0ce001eSMatthew G. Knepley 1096c0ce001eSMatthew G. Knepley Input Parameter: 1097c0ce001eSMatthew G. Knepley . fvm - the PetscFV object to set options for 1098c0ce001eSMatthew G. Knepley 1099c0ce001eSMatthew G. Knepley Options Database Key: 1100c0ce001eSMatthew G. Knepley . -petscfv_compute_gradients <bool> - Determines whether cell gradients are calculated 1101c0ce001eSMatthew G. Knepley 110288f5f89eSMatthew G. Knepley Level: intermediate 1103c0ce001eSMatthew G. Knepley 1104c0ce001eSMatthew G. Knepley .seealso: PetscFVView() 1105c0ce001eSMatthew G. Knepley @*/ 1106c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetFromOptions(PetscFV fvm) 1107c0ce001eSMatthew G. Knepley { 1108c0ce001eSMatthew G. Knepley const char *defaultType; 1109c0ce001eSMatthew G. Knepley char name[256]; 1110c0ce001eSMatthew G. Knepley PetscBool flg; 1111c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1112c0ce001eSMatthew G. Knepley 1113c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1114c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1115c0ce001eSMatthew G. Knepley if (!((PetscObject) fvm)->type_name) defaultType = PETSCFVUPWIND; 1116c0ce001eSMatthew G. Knepley else defaultType = ((PetscObject) fvm)->type_name; 1117c0ce001eSMatthew G. Knepley ierr = PetscFVRegisterAll();CHKERRQ(ierr); 1118c0ce001eSMatthew G. Knepley 1119c0ce001eSMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) fvm);CHKERRQ(ierr); 1120c0ce001eSMatthew G. Knepley ierr = PetscOptionsFList("-petscfv_type", "Finite volume discretization", "PetscFVSetType", PetscFVList, defaultType, name, 256, &flg);CHKERRQ(ierr); 1121c0ce001eSMatthew G. Knepley if (flg) { 1122c0ce001eSMatthew G. Knepley ierr = PetscFVSetType(fvm, name);CHKERRQ(ierr); 1123c0ce001eSMatthew G. Knepley } else if (!((PetscObject) fvm)->type_name) { 1124c0ce001eSMatthew G. Knepley ierr = PetscFVSetType(fvm, defaultType);CHKERRQ(ierr); 1125c0ce001eSMatthew G. Knepley 1126c0ce001eSMatthew G. Knepley } 1127c0ce001eSMatthew G. Knepley ierr = PetscOptionsBool("-petscfv_compute_gradients", "Compute cell gradients", "PetscFVSetComputeGradients", fvm->computeGradients, &fvm->computeGradients, NULL);CHKERRQ(ierr); 1128c0ce001eSMatthew G. Knepley if (fvm->ops->setfromoptions) {ierr = (*fvm->ops->setfromoptions)(fvm);CHKERRQ(ierr);} 1129c0ce001eSMatthew G. Knepley /* process any options handlers added with PetscObjectAddOptionsHandler() */ 1130c0ce001eSMatthew G. Knepley ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fvm);CHKERRQ(ierr); 1131c0ce001eSMatthew G. Knepley ierr = PetscLimiterSetFromOptions(fvm->limiter);CHKERRQ(ierr); 1132c0ce001eSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 1133c0ce001eSMatthew G. Knepley ierr = PetscFVViewFromOptions(fvm, NULL, "-petscfv_view");CHKERRQ(ierr); 1134c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1135c0ce001eSMatthew G. Knepley } 1136c0ce001eSMatthew G. Knepley 1137c0ce001eSMatthew G. Knepley /*@ 1138c0ce001eSMatthew G. Knepley PetscFVSetUp - Construct data structures for the PetscFV 1139c0ce001eSMatthew G. Knepley 1140c0ce001eSMatthew G. Knepley Collective on fvm 1141c0ce001eSMatthew G. Knepley 1142c0ce001eSMatthew G. Knepley Input Parameter: 1143c0ce001eSMatthew G. Knepley . fvm - the PetscFV object to setup 1144c0ce001eSMatthew G. Knepley 114588f5f89eSMatthew G. Knepley Level: intermediate 1146c0ce001eSMatthew G. Knepley 1147c0ce001eSMatthew G. Knepley .seealso: PetscFVView(), PetscFVDestroy() 1148c0ce001eSMatthew G. Knepley @*/ 1149c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetUp(PetscFV fvm) 1150c0ce001eSMatthew G. Knepley { 1151c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1152c0ce001eSMatthew G. Knepley 1153c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1154c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1155c0ce001eSMatthew G. Knepley ierr = PetscLimiterSetUp(fvm->limiter);CHKERRQ(ierr); 1156c0ce001eSMatthew G. Knepley if (fvm->ops->setup) {ierr = (*fvm->ops->setup)(fvm);CHKERRQ(ierr);} 1157c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1158c0ce001eSMatthew G. Knepley } 1159c0ce001eSMatthew G. Knepley 1160c0ce001eSMatthew G. Knepley /*@ 1161c0ce001eSMatthew G. Knepley PetscFVDestroy - Destroys a PetscFV object 1162c0ce001eSMatthew G. Knepley 1163c0ce001eSMatthew G. Knepley Collective on fvm 1164c0ce001eSMatthew G. Knepley 1165c0ce001eSMatthew G. Knepley Input Parameter: 1166c0ce001eSMatthew G. Knepley . fvm - the PetscFV object to destroy 1167c0ce001eSMatthew G. Knepley 116888f5f89eSMatthew G. Knepley Level: beginner 1169c0ce001eSMatthew G. Knepley 1170c0ce001eSMatthew G. Knepley .seealso: PetscFVView() 1171c0ce001eSMatthew G. Knepley @*/ 1172c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVDestroy(PetscFV *fvm) 1173c0ce001eSMatthew G. Knepley { 1174c0ce001eSMatthew G. Knepley PetscInt i; 1175c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1176c0ce001eSMatthew G. Knepley 1177c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1178c0ce001eSMatthew G. Knepley if (!*fvm) PetscFunctionReturn(0); 1179c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific((*fvm), PETSCFV_CLASSID, 1); 1180c0ce001eSMatthew G. Knepley 1181c0ce001eSMatthew G. Knepley if (--((PetscObject)(*fvm))->refct > 0) {*fvm = 0; PetscFunctionReturn(0);} 1182c0ce001eSMatthew G. Knepley ((PetscObject) (*fvm))->refct = 0; 1183c0ce001eSMatthew G. Knepley 1184c0ce001eSMatthew G. Knepley for (i = 0; i < (*fvm)->numComponents; i++) { 1185c0ce001eSMatthew G. Knepley ierr = PetscFree((*fvm)->componentNames[i]);CHKERRQ(ierr); 1186c0ce001eSMatthew G. Knepley } 1187c0ce001eSMatthew G. Knepley ierr = PetscFree((*fvm)->componentNames);CHKERRQ(ierr); 1188c0ce001eSMatthew G. Knepley ierr = PetscLimiterDestroy(&(*fvm)->limiter);CHKERRQ(ierr); 1189c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceDestroy(&(*fvm)->dualSpace);CHKERRQ(ierr); 1190c0ce001eSMatthew G. Knepley ierr = PetscFree((*fvm)->fluxWork);CHKERRQ(ierr); 1191c0ce001eSMatthew G. Knepley ierr = PetscQuadratureDestroy(&(*fvm)->quadrature);CHKERRQ(ierr); 1192ef0bb6c7SMatthew G. Knepley ierr = PetscTabulationDestroy(&(*fvm)->T);CHKERRQ(ierr); 1193c0ce001eSMatthew G. Knepley 1194c0ce001eSMatthew G. Knepley if ((*fvm)->ops->destroy) {ierr = (*(*fvm)->ops->destroy)(*fvm);CHKERRQ(ierr);} 1195c0ce001eSMatthew G. Knepley ierr = PetscHeaderDestroy(fvm);CHKERRQ(ierr); 1196c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1197c0ce001eSMatthew G. Knepley } 1198c0ce001eSMatthew G. Knepley 1199c0ce001eSMatthew G. Knepley /*@ 1200c0ce001eSMatthew G. Knepley PetscFVCreate - Creates an empty PetscFV object. The type can then be set with PetscFVSetType(). 1201c0ce001eSMatthew G. Knepley 1202c0ce001eSMatthew G. Knepley Collective 1203c0ce001eSMatthew G. Knepley 1204c0ce001eSMatthew G. Knepley Input Parameter: 1205c0ce001eSMatthew G. Knepley . comm - The communicator for the PetscFV object 1206c0ce001eSMatthew G. Knepley 1207c0ce001eSMatthew G. Knepley Output Parameter: 1208c0ce001eSMatthew G. Knepley . fvm - The PetscFV object 1209c0ce001eSMatthew G. Knepley 1210c0ce001eSMatthew G. Knepley Level: beginner 1211c0ce001eSMatthew G. Knepley 1212c0ce001eSMatthew G. Knepley .seealso: PetscFVSetType(), PETSCFVUPWIND 1213c0ce001eSMatthew G. Knepley @*/ 1214c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVCreate(MPI_Comm comm, PetscFV *fvm) 1215c0ce001eSMatthew G. Knepley { 1216c0ce001eSMatthew G. Knepley PetscFV f; 1217c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1218c0ce001eSMatthew G. Knepley 1219c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1220c0ce001eSMatthew G. Knepley PetscValidPointer(fvm, 2); 1221c0ce001eSMatthew G. Knepley *fvm = NULL; 1222c0ce001eSMatthew G. Knepley ierr = PetscFVInitializePackage();CHKERRQ(ierr); 1223c0ce001eSMatthew G. Knepley 1224c0ce001eSMatthew G. Knepley ierr = PetscHeaderCreate(f, PETSCFV_CLASSID, "PetscFV", "Finite Volume", "PetscFV", comm, PetscFVDestroy, PetscFVView);CHKERRQ(ierr); 1225c0ce001eSMatthew G. Knepley ierr = PetscMemzero(f->ops, sizeof(struct _PetscFVOps));CHKERRQ(ierr); 1226c0ce001eSMatthew G. Knepley 1227c0ce001eSMatthew G. Knepley ierr = PetscLimiterCreate(comm, &f->limiter);CHKERRQ(ierr); 1228c0ce001eSMatthew G. Knepley f->numComponents = 1; 1229c0ce001eSMatthew G. Knepley f->dim = 0; 1230c0ce001eSMatthew G. Knepley f->computeGradients = PETSC_FALSE; 1231c0ce001eSMatthew G. Knepley f->fluxWork = NULL; 1232c0ce001eSMatthew G. Knepley ierr = PetscCalloc1(f->numComponents,&f->componentNames);CHKERRQ(ierr); 1233c0ce001eSMatthew G. Knepley 1234c0ce001eSMatthew G. Knepley *fvm = f; 1235c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1236c0ce001eSMatthew G. Knepley } 1237c0ce001eSMatthew G. Knepley 1238c0ce001eSMatthew G. Knepley /*@ 1239c0ce001eSMatthew G. Knepley PetscFVSetLimiter - Set the limiter object 1240c0ce001eSMatthew G. Knepley 1241c0ce001eSMatthew G. Knepley Logically collective on fvm 1242c0ce001eSMatthew G. Knepley 1243c0ce001eSMatthew G. Knepley Input Parameters: 1244c0ce001eSMatthew G. Knepley + fvm - the PetscFV object 1245c0ce001eSMatthew G. Knepley - lim - The PetscLimiter 1246c0ce001eSMatthew G. Knepley 124788f5f89eSMatthew G. Knepley Level: intermediate 1248c0ce001eSMatthew G. Knepley 1249c0ce001eSMatthew G. Knepley .seealso: PetscFVGetLimiter() 1250c0ce001eSMatthew G. Knepley @*/ 1251c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetLimiter(PetscFV fvm, PetscLimiter lim) 1252c0ce001eSMatthew G. Knepley { 1253c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1254c0ce001eSMatthew G. Knepley 1255c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1256c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1257c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 2); 1258c0ce001eSMatthew G. Knepley ierr = PetscLimiterDestroy(&fvm->limiter);CHKERRQ(ierr); 1259c0ce001eSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) lim);CHKERRQ(ierr); 1260c0ce001eSMatthew G. Knepley fvm->limiter = lim; 1261c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1262c0ce001eSMatthew G. Knepley } 1263c0ce001eSMatthew G. Knepley 1264c0ce001eSMatthew G. Knepley /*@ 1265c0ce001eSMatthew G. Knepley PetscFVGetLimiter - Get the limiter object 1266c0ce001eSMatthew G. Knepley 1267c0ce001eSMatthew G. Knepley Not collective 1268c0ce001eSMatthew G. Knepley 1269c0ce001eSMatthew G. Knepley Input Parameter: 1270c0ce001eSMatthew G. Knepley . fvm - the PetscFV object 1271c0ce001eSMatthew G. Knepley 1272c0ce001eSMatthew G. Knepley Output Parameter: 1273c0ce001eSMatthew G. Knepley . lim - The PetscLimiter 1274c0ce001eSMatthew G. Knepley 127588f5f89eSMatthew G. Knepley Level: intermediate 1276c0ce001eSMatthew G. Knepley 1277c0ce001eSMatthew G. Knepley .seealso: PetscFVSetLimiter() 1278c0ce001eSMatthew G. Knepley @*/ 1279c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetLimiter(PetscFV fvm, PetscLimiter *lim) 1280c0ce001eSMatthew G. Knepley { 1281c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1282c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1283c0ce001eSMatthew G. Knepley PetscValidPointer(lim, 2); 1284c0ce001eSMatthew G. Knepley *lim = fvm->limiter; 1285c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1286c0ce001eSMatthew G. Knepley } 1287c0ce001eSMatthew G. Knepley 1288c0ce001eSMatthew G. Knepley /*@ 1289c0ce001eSMatthew G. Knepley PetscFVSetNumComponents - Set the number of field components 1290c0ce001eSMatthew G. Knepley 1291c0ce001eSMatthew G. Knepley Logically collective on fvm 1292c0ce001eSMatthew G. Knepley 1293c0ce001eSMatthew G. Knepley Input Parameters: 1294c0ce001eSMatthew G. Knepley + fvm - the PetscFV object 1295c0ce001eSMatthew G. Knepley - comp - The number of components 1296c0ce001eSMatthew G. Knepley 129788f5f89eSMatthew G. Knepley Level: intermediate 1298c0ce001eSMatthew G. Knepley 1299c0ce001eSMatthew G. Knepley .seealso: PetscFVGetNumComponents() 1300c0ce001eSMatthew G. Knepley @*/ 1301c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetNumComponents(PetscFV fvm, PetscInt comp) 1302c0ce001eSMatthew G. Knepley { 1303c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1304c0ce001eSMatthew G. Knepley 1305c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1306c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1307c0ce001eSMatthew G. Knepley if (fvm->numComponents != comp) { 1308c0ce001eSMatthew G. Knepley PetscInt i; 1309c0ce001eSMatthew G. Knepley 1310c0ce001eSMatthew G. Knepley for (i = 0; i < fvm->numComponents; i++) { 1311c0ce001eSMatthew G. Knepley ierr = PetscFree(fvm->componentNames[i]);CHKERRQ(ierr); 1312c0ce001eSMatthew G. Knepley } 1313c0ce001eSMatthew G. Knepley ierr = PetscFree(fvm->componentNames);CHKERRQ(ierr); 1314c0ce001eSMatthew G. Knepley ierr = PetscCalloc1(comp,&fvm->componentNames);CHKERRQ(ierr); 1315c0ce001eSMatthew G. Knepley } 1316c0ce001eSMatthew G. Knepley fvm->numComponents = comp; 1317c0ce001eSMatthew G. Knepley ierr = PetscFree(fvm->fluxWork);CHKERRQ(ierr); 1318c0ce001eSMatthew G. Knepley ierr = PetscMalloc1(comp, &fvm->fluxWork);CHKERRQ(ierr); 1319c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1320c0ce001eSMatthew G. Knepley } 1321c0ce001eSMatthew G. Knepley 1322c0ce001eSMatthew G. Knepley /*@ 1323c0ce001eSMatthew G. Knepley PetscFVGetNumComponents - Get the number of field components 1324c0ce001eSMatthew G. Knepley 1325c0ce001eSMatthew G. Knepley Not collective 1326c0ce001eSMatthew G. Knepley 1327c0ce001eSMatthew G. Knepley Input Parameter: 1328c0ce001eSMatthew G. Knepley . fvm - the PetscFV object 1329c0ce001eSMatthew G. Knepley 1330c0ce001eSMatthew G. Knepley Output Parameter: 1331c0ce001eSMatthew G. Knepley , comp - The number of components 1332c0ce001eSMatthew G. Knepley 133388f5f89eSMatthew G. Knepley Level: intermediate 1334c0ce001eSMatthew G. Knepley 1335c0ce001eSMatthew G. Knepley .seealso: PetscFVSetNumComponents() 1336c0ce001eSMatthew G. Knepley @*/ 1337c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetNumComponents(PetscFV fvm, PetscInt *comp) 1338c0ce001eSMatthew G. Knepley { 1339c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1340c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1341c0ce001eSMatthew G. Knepley PetscValidPointer(comp, 2); 1342c0ce001eSMatthew G. Knepley *comp = fvm->numComponents; 1343c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1344c0ce001eSMatthew G. Knepley } 1345c0ce001eSMatthew G. Knepley 1346c0ce001eSMatthew G. Knepley /*@C 1347c0ce001eSMatthew G. Knepley PetscFVSetComponentName - Set the name of a component (used in output and viewing) 1348c0ce001eSMatthew G. Knepley 1349c0ce001eSMatthew G. Knepley Logically collective on fvm 1350c0ce001eSMatthew G. Knepley Input Parameters: 1351c0ce001eSMatthew G. Knepley + fvm - the PetscFV object 1352c0ce001eSMatthew G. Knepley . comp - the component number 1353c0ce001eSMatthew G. Knepley - name - the component name 1354c0ce001eSMatthew G. Knepley 135588f5f89eSMatthew G. Knepley Level: intermediate 1356c0ce001eSMatthew G. Knepley 1357c0ce001eSMatthew G. Knepley .seealso: PetscFVGetComponentName() 1358c0ce001eSMatthew G. Knepley @*/ 1359c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetComponentName(PetscFV fvm, PetscInt comp, const char *name) 1360c0ce001eSMatthew G. Knepley { 1361c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1362c0ce001eSMatthew G. Knepley 1363c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1364c0ce001eSMatthew G. Knepley ierr = PetscFree(fvm->componentNames[comp]);CHKERRQ(ierr); 1365c0ce001eSMatthew G. Knepley ierr = PetscStrallocpy(name,&fvm->componentNames[comp]);CHKERRQ(ierr); 1366c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1367c0ce001eSMatthew G. Knepley } 1368c0ce001eSMatthew G. Knepley 1369c0ce001eSMatthew G. Knepley /*@C 1370c0ce001eSMatthew G. Knepley PetscFVGetComponentName - Get the name of a component (used in output and viewing) 1371c0ce001eSMatthew G. Knepley 1372c0ce001eSMatthew G. Knepley Logically collective on fvm 1373c0ce001eSMatthew G. Knepley Input Parameters: 1374c0ce001eSMatthew G. Knepley + fvm - the PetscFV object 1375c0ce001eSMatthew G. Knepley - comp - the component number 1376c0ce001eSMatthew G. Knepley 1377c0ce001eSMatthew G. Knepley Output Parameter: 1378c0ce001eSMatthew G. Knepley . name - the component name 1379c0ce001eSMatthew G. Knepley 138088f5f89eSMatthew G. Knepley Level: intermediate 1381c0ce001eSMatthew G. Knepley 1382c0ce001eSMatthew G. Knepley .seealso: PetscFVSetComponentName() 1383c0ce001eSMatthew G. Knepley @*/ 1384c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetComponentName(PetscFV fvm, PetscInt comp, const char **name) 1385c0ce001eSMatthew G. Knepley { 1386c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1387c0ce001eSMatthew G. Knepley *name = fvm->componentNames[comp]; 1388c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1389c0ce001eSMatthew G. Knepley } 1390c0ce001eSMatthew G. Knepley 1391c0ce001eSMatthew G. Knepley /*@ 1392c0ce001eSMatthew G. Knepley PetscFVSetSpatialDimension - Set the spatial dimension 1393c0ce001eSMatthew G. Knepley 1394c0ce001eSMatthew G. Knepley Logically collective on fvm 1395c0ce001eSMatthew G. Knepley 1396c0ce001eSMatthew G. Knepley Input Parameters: 1397c0ce001eSMatthew G. Knepley + fvm - the PetscFV object 1398c0ce001eSMatthew G. Knepley - dim - The spatial dimension 1399c0ce001eSMatthew G. Knepley 140088f5f89eSMatthew G. Knepley Level: intermediate 1401c0ce001eSMatthew G. Knepley 1402c0ce001eSMatthew G. Knepley .seealso: PetscFVGetSpatialDimension() 1403c0ce001eSMatthew G. Knepley @*/ 1404c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetSpatialDimension(PetscFV fvm, PetscInt dim) 1405c0ce001eSMatthew G. Knepley { 1406c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1407c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1408c0ce001eSMatthew G. Knepley fvm->dim = dim; 1409c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1410c0ce001eSMatthew G. Knepley } 1411c0ce001eSMatthew G. Knepley 1412c0ce001eSMatthew G. Knepley /*@ 1413c0ce001eSMatthew G. Knepley PetscFVGetSpatialDimension - Get the spatial dimension 1414c0ce001eSMatthew G. Knepley 1415c0ce001eSMatthew G. Knepley Logically collective on fvm 1416c0ce001eSMatthew G. Knepley 1417c0ce001eSMatthew G. Knepley Input Parameter: 1418c0ce001eSMatthew G. Knepley . fvm - the PetscFV object 1419c0ce001eSMatthew G. Knepley 1420c0ce001eSMatthew G. Knepley Output Parameter: 1421c0ce001eSMatthew G. Knepley . dim - The spatial dimension 1422c0ce001eSMatthew G. Knepley 142388f5f89eSMatthew G. Knepley Level: intermediate 1424c0ce001eSMatthew G. Knepley 1425c0ce001eSMatthew G. Knepley .seealso: PetscFVSetSpatialDimension() 1426c0ce001eSMatthew G. Knepley @*/ 1427c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetSpatialDimension(PetscFV fvm, PetscInt *dim) 1428c0ce001eSMatthew G. Knepley { 1429c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1430c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1431c0ce001eSMatthew G. Knepley PetscValidPointer(dim, 2); 1432c0ce001eSMatthew G. Knepley *dim = fvm->dim; 1433c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1434c0ce001eSMatthew G. Knepley } 1435c0ce001eSMatthew G. Knepley 1436c0ce001eSMatthew G. Knepley /*@ 1437c0ce001eSMatthew G. Knepley PetscFVSetComputeGradients - Toggle computation of cell gradients 1438c0ce001eSMatthew G. Knepley 1439c0ce001eSMatthew G. Knepley Logically collective on fvm 1440c0ce001eSMatthew G. Knepley 1441c0ce001eSMatthew G. Knepley Input Parameters: 1442c0ce001eSMatthew G. Knepley + fvm - the PetscFV object 1443c0ce001eSMatthew G. Knepley - computeGradients - Flag to compute cell gradients 1444c0ce001eSMatthew G. Knepley 144588f5f89eSMatthew G. Knepley Level: intermediate 1446c0ce001eSMatthew G. Knepley 1447c0ce001eSMatthew G. Knepley .seealso: PetscFVGetComputeGradients() 1448c0ce001eSMatthew G. Knepley @*/ 1449c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetComputeGradients(PetscFV fvm, PetscBool computeGradients) 1450c0ce001eSMatthew G. Knepley { 1451c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1452c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1453c0ce001eSMatthew G. Knepley fvm->computeGradients = computeGradients; 1454c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1455c0ce001eSMatthew G. Knepley } 1456c0ce001eSMatthew G. Knepley 1457c0ce001eSMatthew G. Knepley /*@ 1458c0ce001eSMatthew G. Knepley PetscFVGetComputeGradients - Return flag for computation of cell gradients 1459c0ce001eSMatthew G. Knepley 1460c0ce001eSMatthew G. Knepley Not collective 1461c0ce001eSMatthew G. Knepley 1462c0ce001eSMatthew G. Knepley Input Parameter: 1463c0ce001eSMatthew G. Knepley . fvm - the PetscFV object 1464c0ce001eSMatthew G. Knepley 1465c0ce001eSMatthew G. Knepley Output Parameter: 1466c0ce001eSMatthew G. Knepley . computeGradients - Flag to compute cell gradients 1467c0ce001eSMatthew G. Knepley 146888f5f89eSMatthew G. Knepley Level: intermediate 1469c0ce001eSMatthew G. Knepley 1470c0ce001eSMatthew G. Knepley .seealso: PetscFVSetComputeGradients() 1471c0ce001eSMatthew G. Knepley @*/ 1472c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetComputeGradients(PetscFV fvm, PetscBool *computeGradients) 1473c0ce001eSMatthew G. Knepley { 1474c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1475c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1476c0ce001eSMatthew G. Knepley PetscValidPointer(computeGradients, 2); 1477c0ce001eSMatthew G. Knepley *computeGradients = fvm->computeGradients; 1478c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1479c0ce001eSMatthew G. Knepley } 1480c0ce001eSMatthew G. Knepley 1481c0ce001eSMatthew G. Knepley /*@ 1482c0ce001eSMatthew G. Knepley PetscFVSetQuadrature - Set the quadrature object 1483c0ce001eSMatthew G. Knepley 1484c0ce001eSMatthew G. Knepley Logically collective on fvm 1485c0ce001eSMatthew G. Knepley 1486c0ce001eSMatthew G. Knepley Input Parameters: 1487c0ce001eSMatthew G. Knepley + fvm - the PetscFV object 1488c0ce001eSMatthew G. Knepley - q - The PetscQuadrature 1489c0ce001eSMatthew G. Knepley 149088f5f89eSMatthew G. Knepley Level: intermediate 1491c0ce001eSMatthew G. Knepley 1492c0ce001eSMatthew G. Knepley .seealso: PetscFVGetQuadrature() 1493c0ce001eSMatthew G. Knepley @*/ 1494c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetQuadrature(PetscFV fvm, PetscQuadrature q) 1495c0ce001eSMatthew G. Knepley { 1496c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1497c0ce001eSMatthew G. Knepley 1498c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1499c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1500c0ce001eSMatthew G. Knepley ierr = PetscQuadratureDestroy(&fvm->quadrature);CHKERRQ(ierr); 1501c0ce001eSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr); 1502c0ce001eSMatthew G. Knepley fvm->quadrature = q; 1503c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1504c0ce001eSMatthew G. Knepley } 1505c0ce001eSMatthew G. Knepley 1506c0ce001eSMatthew G. Knepley /*@ 1507c0ce001eSMatthew G. Knepley PetscFVGetQuadrature - Get the quadrature object 1508c0ce001eSMatthew G. Knepley 1509c0ce001eSMatthew G. Knepley Not collective 1510c0ce001eSMatthew G. Knepley 1511c0ce001eSMatthew G. Knepley Input Parameter: 1512c0ce001eSMatthew G. Knepley . fvm - the PetscFV object 1513c0ce001eSMatthew G. Knepley 1514c0ce001eSMatthew G. Knepley Output Parameter: 1515c0ce001eSMatthew G. Knepley . lim - The PetscQuadrature 1516c0ce001eSMatthew G. Knepley 151788f5f89eSMatthew G. Knepley Level: intermediate 1518c0ce001eSMatthew G. Knepley 1519c0ce001eSMatthew G. Knepley .seealso: PetscFVSetQuadrature() 1520c0ce001eSMatthew G. Knepley @*/ 1521c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetQuadrature(PetscFV fvm, PetscQuadrature *q) 1522c0ce001eSMatthew G. Knepley { 1523c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1524c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1525c0ce001eSMatthew G. Knepley PetscValidPointer(q, 2); 1526c0ce001eSMatthew G. Knepley if (!fvm->quadrature) { 1527c0ce001eSMatthew G. Knepley /* Create default 1-point quadrature */ 1528c0ce001eSMatthew G. Knepley PetscReal *points, *weights; 1529c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1530c0ce001eSMatthew G. Knepley 1531c0ce001eSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &fvm->quadrature);CHKERRQ(ierr); 1532c0ce001eSMatthew G. Knepley ierr = PetscCalloc1(fvm->dim, &points);CHKERRQ(ierr); 1533c0ce001eSMatthew G. Knepley ierr = PetscMalloc1(1, &weights);CHKERRQ(ierr); 1534c0ce001eSMatthew G. Knepley weights[0] = 1.0; 1535c0ce001eSMatthew G. Knepley ierr = PetscQuadratureSetData(fvm->quadrature, fvm->dim, 1, 1, points, weights);CHKERRQ(ierr); 1536c0ce001eSMatthew G. Knepley } 1537c0ce001eSMatthew G. Knepley *q = fvm->quadrature; 1538c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1539c0ce001eSMatthew G. Knepley } 1540c0ce001eSMatthew G. Knepley 1541c0ce001eSMatthew G. Knepley /*@ 1542c0ce001eSMatthew G. Knepley PetscFVGetDualSpace - Returns the PetscDualSpace used to define the inner product 1543c0ce001eSMatthew G. Knepley 1544c0ce001eSMatthew G. Knepley Not collective 1545c0ce001eSMatthew G. Knepley 1546c0ce001eSMatthew G. Knepley Input Parameter: 1547c0ce001eSMatthew G. Knepley . fvm - The PetscFV object 1548c0ce001eSMatthew G. Knepley 1549c0ce001eSMatthew G. Knepley Output Parameter: 1550c0ce001eSMatthew G. Knepley . sp - The PetscDualSpace object 1551c0ce001eSMatthew G. Knepley 1552c0ce001eSMatthew G. Knepley Note: A simple dual space is provided automatically, and the user typically will not need to override it. 1553c0ce001eSMatthew G. Knepley 155488f5f89eSMatthew G. Knepley Level: intermediate 1555c0ce001eSMatthew G. Knepley 1556c0ce001eSMatthew G. Knepley .seealso: PetscFVCreate() 1557c0ce001eSMatthew G. Knepley @*/ 1558c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetDualSpace(PetscFV fvm, PetscDualSpace *sp) 1559c0ce001eSMatthew G. Knepley { 1560c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1561c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1562c0ce001eSMatthew G. Knepley PetscValidPointer(sp, 2); 1563c0ce001eSMatthew G. Knepley if (!fvm->dualSpace) { 1564c0ce001eSMatthew G. Knepley DM K; 1565c0ce001eSMatthew G. Knepley PetscInt dim, Nc, c; 1566c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1567c0ce001eSMatthew G. Knepley 1568c0ce001eSMatthew G. Knepley ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr); 1569c0ce001eSMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &Nc);CHKERRQ(ierr); 1570c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fvm), &fvm->dualSpace);CHKERRQ(ierr); 1571c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceSetType(fvm->dualSpace, PETSCDUALSPACESIMPLE);CHKERRQ(ierr); 1572c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceCreateReferenceCell(fvm->dualSpace, dim, PETSC_FALSE, &K);CHKERRQ(ierr); /* TODO: The reference cell type should be held by the discretization object */ 1573c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceSetNumComponents(fvm->dualSpace, Nc);CHKERRQ(ierr); 1574c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceSetDM(fvm->dualSpace, K);CHKERRQ(ierr); 1575c0ce001eSMatthew G. Knepley ierr = DMDestroy(&K);CHKERRQ(ierr); 1576c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceSimpleSetDimension(fvm->dualSpace, Nc);CHKERRQ(ierr); 1577c0ce001eSMatthew G. Knepley /* Should we be using PetscFVGetQuadrature() here? */ 1578c0ce001eSMatthew G. Knepley for (c = 0; c < Nc; ++c) { 1579c0ce001eSMatthew G. Knepley PetscQuadrature qc; 1580c0ce001eSMatthew G. Knepley PetscReal *points, *weights; 1581c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1582c0ce001eSMatthew G. Knepley 1583c0ce001eSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &qc);CHKERRQ(ierr); 1584c0ce001eSMatthew G. Knepley ierr = PetscCalloc1(dim, &points);CHKERRQ(ierr); 1585c0ce001eSMatthew G. Knepley ierr = PetscCalloc1(Nc, &weights);CHKERRQ(ierr); 1586c0ce001eSMatthew G. Knepley weights[c] = 1.0; 1587c0ce001eSMatthew G. Knepley ierr = PetscQuadratureSetData(qc, dim, Nc, 1, points, weights);CHKERRQ(ierr); 1588c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceSimpleSetFunctional(fvm->dualSpace, c, qc);CHKERRQ(ierr); 1589c0ce001eSMatthew G. Knepley ierr = PetscQuadratureDestroy(&qc);CHKERRQ(ierr); 1590c0ce001eSMatthew G. Knepley } 1591c0ce001eSMatthew G. Knepley } 1592c0ce001eSMatthew G. Knepley *sp = fvm->dualSpace; 1593c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1594c0ce001eSMatthew G. Knepley } 1595c0ce001eSMatthew G. Knepley 1596c0ce001eSMatthew G. Knepley /*@ 1597c0ce001eSMatthew G. Knepley PetscFVSetDualSpace - Sets the PetscDualSpace used to define the inner product 1598c0ce001eSMatthew G. Knepley 1599c0ce001eSMatthew G. Knepley Not collective 1600c0ce001eSMatthew G. Knepley 1601c0ce001eSMatthew G. Knepley Input Parameters: 1602c0ce001eSMatthew G. Knepley + fvm - The PetscFV object 1603c0ce001eSMatthew G. Knepley - sp - The PetscDualSpace object 1604c0ce001eSMatthew G. Knepley 1605c0ce001eSMatthew G. Knepley Level: intermediate 1606c0ce001eSMatthew G. Knepley 1607c0ce001eSMatthew G. Knepley Note: A simple dual space is provided automatically, and the user typically will not need to override it. 1608c0ce001eSMatthew G. Knepley 1609c0ce001eSMatthew G. Knepley .seealso: PetscFVCreate() 1610c0ce001eSMatthew G. Knepley @*/ 1611c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetDualSpace(PetscFV fvm, PetscDualSpace sp) 1612c0ce001eSMatthew G. Knepley { 1613c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1614c0ce001eSMatthew G. Knepley 1615c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1616c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1617c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2); 1618c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceDestroy(&fvm->dualSpace);CHKERRQ(ierr); 1619c0ce001eSMatthew G. Knepley fvm->dualSpace = sp; 1620c0ce001eSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) fvm->dualSpace);CHKERRQ(ierr); 1621c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1622c0ce001eSMatthew G. Knepley } 1623c0ce001eSMatthew G. Knepley 162488f5f89eSMatthew G. Knepley /*@C 1625ef0bb6c7SMatthew G. Knepley PetscFVGetCellTabulation - Returns the tabulation of the basis functions at the quadrature points 162688f5f89eSMatthew G. Knepley 162788f5f89eSMatthew G. Knepley Not collective 162888f5f89eSMatthew G. Knepley 162988f5f89eSMatthew G. Knepley Input Parameter: 163088f5f89eSMatthew G. Knepley . fvm - The PetscFV object 163188f5f89eSMatthew G. Knepley 1632ef0bb6c7SMatthew G. Knepley Output Parameter: 1633ef0bb6c7SMatthew G. Knepley . T - The basis function values and derviatives at quadrature points 163488f5f89eSMatthew G. Knepley 163588f5f89eSMatthew G. Knepley Note: 1636ef0bb6c7SMatthew G. Knepley $ T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 1637ef0bb6c7SMatthew G. Knepley $ T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d 1638ef0bb6c7SMatthew G. Knepley $ T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e 163988f5f89eSMatthew G. Knepley 164088f5f89eSMatthew G. Knepley Level: intermediate 164188f5f89eSMatthew G. Knepley 1642ef0bb6c7SMatthew G. Knepley .seealso: PetscFEGetCellTabulation(), PetscFVCreateTabulation(), PetscFVGetQuadrature(), PetscQuadratureGetData() 164388f5f89eSMatthew G. Knepley @*/ 1644ef0bb6c7SMatthew G. Knepley PetscErrorCode PetscFVGetCellTabulation(PetscFV fvm, PetscTabulation *T) 1645c0ce001eSMatthew G. Knepley { 1646c0ce001eSMatthew G. Knepley PetscInt npoints; 1647c0ce001eSMatthew G. Knepley const PetscReal *points; 1648c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1649c0ce001eSMatthew G. Knepley 1650c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1651c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1652ef0bb6c7SMatthew G. Knepley PetscValidPointer(T, 2); 1653c0ce001eSMatthew G. Knepley ierr = PetscQuadratureGetData(fvm->quadrature, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr); 1654ef0bb6c7SMatthew G. Knepley if (!fvm->T) {ierr = PetscFVCreateTabulation(fvm, 1, npoints, points, 1, &fvm->T);CHKERRQ(ierr);} 1655ef0bb6c7SMatthew G. Knepley *T = fvm->T; 1656c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1657c0ce001eSMatthew G. Knepley } 1658c0ce001eSMatthew G. Knepley 165988f5f89eSMatthew G. Knepley /*@C 1660ef0bb6c7SMatthew G. Knepley PetscFVCreateTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 166188f5f89eSMatthew G. Knepley 166288f5f89eSMatthew G. Knepley Not collective 166388f5f89eSMatthew G. Knepley 166488f5f89eSMatthew G. Knepley Input Parameters: 166588f5f89eSMatthew G. Knepley + fvm - The PetscFV object 1666ef0bb6c7SMatthew G. Knepley . nrepl - The number of replicas 1667ef0bb6c7SMatthew G. Knepley . npoints - The number of tabulation points in a replica 1668ef0bb6c7SMatthew G. Knepley . points - The tabulation point coordinates 1669ef0bb6c7SMatthew G. Knepley - K - The order of derivative to tabulate 167088f5f89eSMatthew G. Knepley 1671ef0bb6c7SMatthew G. Knepley Output Parameter: 1672ef0bb6c7SMatthew G. Knepley . T - The basis function values and derviative at tabulation points 167388f5f89eSMatthew G. Knepley 167488f5f89eSMatthew G. Knepley Note: 1675ef0bb6c7SMatthew G. Knepley $ T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 1676ef0bb6c7SMatthew G. Knepley $ T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d 1677ef0bb6c7SMatthew G. Knepley $ T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e 167888f5f89eSMatthew G. Knepley 167988f5f89eSMatthew G. Knepley Level: intermediate 168088f5f89eSMatthew G. Knepley 1681ef0bb6c7SMatthew G. Knepley .seealso: PetscFECreateTabulation(), PetscTabulationDestroy(), PetscFEGetCellTabulation() 168288f5f89eSMatthew G. Knepley @*/ 1683ef0bb6c7SMatthew G. Knepley PetscErrorCode PetscFVCreateTabulation(PetscFV fvm, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T) 1684c0ce001eSMatthew G. Knepley { 1685c0ce001eSMatthew G. Knepley PetscInt pdim = 1; /* Dimension of approximation space P */ 1686ef0bb6c7SMatthew G. Knepley PetscInt cdim; /* Spatial dimension */ 1687ef0bb6c7SMatthew G. Knepley PetscInt Nc; /* Field components */ 1688ef0bb6c7SMatthew G. Knepley PetscInt k, p, d, c, e; 1689c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1690c0ce001eSMatthew G. Knepley 1691c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1692ef0bb6c7SMatthew G. Knepley if (!npoints || K < 0) { 1693ef0bb6c7SMatthew G. Knepley *T = NULL; 1694c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1695c0ce001eSMatthew G. Knepley } 1696c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 169740a2aa30SMatthew G. Knepley PetscValidPointer(points, 4); 169840a2aa30SMatthew G. Knepley PetscValidPointer(T, 6); 1699ef0bb6c7SMatthew G. Knepley ierr = PetscFVGetSpatialDimension(fvm, &cdim);CHKERRQ(ierr); 1700ef0bb6c7SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &Nc);CHKERRQ(ierr); 1701ef0bb6c7SMatthew G. Knepley ierr = PetscMalloc1(1, T);CHKERRQ(ierr); 1702ef0bb6c7SMatthew G. Knepley (*T)->K = !cdim ? 0 : K; 1703ef0bb6c7SMatthew G. Knepley (*T)->Nr = nrepl; 1704ef0bb6c7SMatthew G. Knepley (*T)->Np = npoints; 1705ef0bb6c7SMatthew G. Knepley (*T)->Nb = pdim; 1706ef0bb6c7SMatthew G. Knepley (*T)->Nc = Nc; 1707ef0bb6c7SMatthew G. Knepley (*T)->cdim = cdim; 1708ef0bb6c7SMatthew G. Knepley ierr = PetscMalloc1((*T)->K+1, &(*T)->T);CHKERRQ(ierr); 1709ef0bb6c7SMatthew G. Knepley for (k = 0; k <= (*T)->K; ++k) { 1710ef0bb6c7SMatthew G. Knepley ierr = PetscMalloc1(nrepl*npoints*pdim*Nc*PetscPowInt(cdim, k), &(*T)->T[k]);CHKERRQ(ierr); 1711ef0bb6c7SMatthew G. Knepley } 1712ef0bb6c7SMatthew G. Knepley if (K >= 0) {for (p = 0; p < nrepl*npoints; ++p) for (d = 0; d < pdim; ++d) for (c = 0; c < Nc; ++c) (*T)->T[0][(p*pdim + d)*Nc + c] = 1.0;} 1713ef0bb6c7SMatthew G. Knepley if (K >= 1) {for (p = 0; p < nrepl*npoints; ++p) for (d = 0; d < pdim; ++d) for (c = 0; c < Nc; ++c) for (e = 0; e < cdim; ++e) (*T)->T[1][((p*pdim + d)*Nc + c)*cdim + e] = 0.0;} 1714ef0bb6c7SMatthew G. Knepley if (K >= 2) {for (p = 0; p < nrepl*npoints; ++p) for (d = 0; d < pdim; ++d) for (c = 0; c < Nc; ++c) for (e = 0; e < cdim*cdim; ++e) (*T)->T[2][((p*pdim + d)*Nc + c)*cdim*cdim + e] = 0.0;} 1715c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1716c0ce001eSMatthew G. Knepley } 1717c0ce001eSMatthew G. Knepley 1718c0ce001eSMatthew G. Knepley /*@C 1719c0ce001eSMatthew G. Knepley PetscFVComputeGradient - Compute the gradient reconstruction matrix for a given cell 1720c0ce001eSMatthew G. Knepley 1721c0ce001eSMatthew G. Knepley Input Parameters: 1722c0ce001eSMatthew G. Knepley + fvm - The PetscFV object 1723c0ce001eSMatthew G. Knepley . numFaces - The number of cell faces which are not constrained 1724c0ce001eSMatthew G. Knepley - dx - The vector from the cell centroid to the neighboring cell centroid for each face 1725c0ce001eSMatthew G. Knepley 172688f5f89eSMatthew G. Knepley Level: advanced 1727c0ce001eSMatthew G. Knepley 1728c0ce001eSMatthew G. Knepley .seealso: PetscFVCreate() 1729c0ce001eSMatthew G. Knepley @*/ 1730c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVComputeGradient(PetscFV fvm, PetscInt numFaces, PetscScalar dx[], PetscScalar grad[]) 1731c0ce001eSMatthew G. Knepley { 1732c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1733c0ce001eSMatthew G. Knepley 1734c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1735c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1736c0ce001eSMatthew G. Knepley if (fvm->ops->computegradient) {ierr = (*fvm->ops->computegradient)(fvm, numFaces, dx, grad);CHKERRQ(ierr);} 1737c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1738c0ce001eSMatthew G. Knepley } 1739c0ce001eSMatthew G. Knepley 174088f5f89eSMatthew G. Knepley /*@C 1741c0ce001eSMatthew G. Knepley PetscFVIntegrateRHSFunction - Produce the cell residual vector for a chunk of elements by quadrature integration 1742c0ce001eSMatthew G. Knepley 1743c0ce001eSMatthew G. Knepley Not collective 1744c0ce001eSMatthew G. Knepley 1745c0ce001eSMatthew G. Knepley Input Parameters: 1746c0ce001eSMatthew G. Knepley + fvm - The PetscFV object for the field being integrated 1747c0ce001eSMatthew G. Knepley . prob - The PetscDS specifing the discretizations and continuum functions 1748c0ce001eSMatthew G. Knepley . field - The field being integrated 1749c0ce001eSMatthew G. Knepley . Nf - The number of faces in the chunk 1750c0ce001eSMatthew G. Knepley . fgeom - The face geometry for each face in the chunk 1751c0ce001eSMatthew G. Knepley . neighborVol - The volume for each pair of cells in the chunk 1752c0ce001eSMatthew G. Knepley . uL - The state from the cell on the left 1753c0ce001eSMatthew G. Knepley - uR - The state from the cell on the right 1754c0ce001eSMatthew G. Knepley 17557a7aea1fSJed Brown Output Parameter: 1756c0ce001eSMatthew G. Knepley + fluxL - the left fluxes for each face 1757c0ce001eSMatthew G. Knepley - fluxR - the right fluxes for each face 175888f5f89eSMatthew G. Knepley 175988f5f89eSMatthew G. Knepley Level: developer 176088f5f89eSMatthew G. Knepley 176188f5f89eSMatthew G. Knepley .seealso: PetscFVCreate() 176288f5f89eSMatthew G. Knepley @*/ 1763c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVIntegrateRHSFunction(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol, 1764c0ce001eSMatthew G. Knepley PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[]) 1765c0ce001eSMatthew G. Knepley { 1766c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1767c0ce001eSMatthew G. Knepley 1768c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1769c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1770c0ce001eSMatthew G. Knepley if (fvm->ops->integraterhsfunction) {ierr = (*fvm->ops->integraterhsfunction)(fvm, prob, field, Nf, fgeom, neighborVol, uL, uR, fluxL, fluxR);CHKERRQ(ierr);} 1771c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1772c0ce001eSMatthew G. Knepley } 1773c0ce001eSMatthew G. Knepley 1774c0ce001eSMatthew G. Knepley /*@ 1775c0ce001eSMatthew G. Knepley PetscFVRefine - Create a "refined" PetscFV object that refines the reference cell into smaller copies. This is typically used 1776c0ce001eSMatthew G. Knepley to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more 1777c0ce001eSMatthew G. Knepley sparsity). It is also used to create an interpolation between regularly refined meshes. 1778c0ce001eSMatthew G. Knepley 1779c0ce001eSMatthew G. Knepley Input Parameter: 1780c0ce001eSMatthew G. Knepley . fv - The initial PetscFV 1781c0ce001eSMatthew G. Knepley 1782c0ce001eSMatthew G. Knepley Output Parameter: 1783c0ce001eSMatthew G. Knepley . fvRef - The refined PetscFV 1784c0ce001eSMatthew G. Knepley 178588f5f89eSMatthew G. Knepley Level: advanced 1786c0ce001eSMatthew G. Knepley 1787c0ce001eSMatthew G. Knepley .seealso: PetscFVType, PetscFVCreate(), PetscFVSetType() 1788c0ce001eSMatthew G. Knepley @*/ 1789c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVRefine(PetscFV fv, PetscFV *fvRef) 1790c0ce001eSMatthew G. Knepley { 1791c0ce001eSMatthew G. Knepley PetscDualSpace Q, Qref; 1792c0ce001eSMatthew G. Knepley DM K, Kref; 1793c0ce001eSMatthew G. Knepley PetscQuadrature q, qref; 1794*412e9a14SMatthew G. Knepley DMPolytopeType ct; 1795*412e9a14SMatthew G. Knepley DMPlexCellRefiner cr; 1796c0ce001eSMatthew G. Knepley PetscReal *v0; 1797c0ce001eSMatthew G. Knepley PetscReal *jac, *invjac; 1798c0ce001eSMatthew G. Knepley PetscInt numComp, numSubelements, s; 1799c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1800c0ce001eSMatthew G. Knepley 1801c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1802c0ce001eSMatthew G. Knepley ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr); 1803c0ce001eSMatthew G. Knepley ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr); 1804c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceGetDM(Q, &K);CHKERRQ(ierr); 1805c0ce001eSMatthew G. Knepley /* Create dual space */ 1806c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceDuplicate(Q, &Qref);CHKERRQ(ierr); 1807c0ce001eSMatthew G. Knepley ierr = DMRefine(K, PetscObjectComm((PetscObject) fv), &Kref);CHKERRQ(ierr); 1808c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceSetDM(Qref, Kref);CHKERRQ(ierr); 1809c0ce001eSMatthew G. Knepley ierr = DMDestroy(&Kref);CHKERRQ(ierr); 1810c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceSetUp(Qref);CHKERRQ(ierr); 1811c0ce001eSMatthew G. Knepley /* Create volume */ 1812c0ce001eSMatthew G. Knepley ierr = PetscFVCreate(PetscObjectComm((PetscObject) fv), fvRef);CHKERRQ(ierr); 1813c0ce001eSMatthew G. Knepley ierr = PetscFVSetDualSpace(*fvRef, Qref);CHKERRQ(ierr); 1814c0ce001eSMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &numComp);CHKERRQ(ierr); 1815c0ce001eSMatthew G. Knepley ierr = PetscFVSetNumComponents(*fvRef, numComp);CHKERRQ(ierr); 1816c0ce001eSMatthew G. Knepley ierr = PetscFVSetUp(*fvRef);CHKERRQ(ierr); 1817c0ce001eSMatthew G. Knepley /* Create quadrature */ 1818*412e9a14SMatthew G. Knepley ierr = DMPlexGetCellType(K, 0, &ct);CHKERRQ(ierr); 1819*412e9a14SMatthew G. Knepley ierr = DMPlexCellRefinerCreate(K, &cr);CHKERRQ(ierr); 1820*412e9a14SMatthew G. Knepley ierr = DMPlexCellRefinerGetAffineTransforms(cr, ct, &numSubelements, &v0, &jac, &invjac);CHKERRQ(ierr); 1821c0ce001eSMatthew G. Knepley ierr = PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);CHKERRQ(ierr); 1822c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceSimpleSetDimension(Qref, numSubelements);CHKERRQ(ierr); 1823c0ce001eSMatthew G. Knepley for (s = 0; s < numSubelements; ++s) { 1824c0ce001eSMatthew G. Knepley PetscQuadrature qs; 1825c0ce001eSMatthew G. Knepley const PetscReal *points, *weights; 1826c0ce001eSMatthew G. Knepley PetscReal *p, *w; 1827c0ce001eSMatthew G. Knepley PetscInt dim, Nc, npoints, np; 1828c0ce001eSMatthew G. Knepley 1829c0ce001eSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &qs);CHKERRQ(ierr); 1830c0ce001eSMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr); 1831c0ce001eSMatthew G. Knepley np = npoints/numSubelements; 1832c0ce001eSMatthew G. Knepley ierr = PetscMalloc1(np*dim,&p);CHKERRQ(ierr); 1833c0ce001eSMatthew G. Knepley ierr = PetscMalloc1(np*Nc,&w);CHKERRQ(ierr); 1834c0ce001eSMatthew G. Knepley ierr = PetscArraycpy(p, &points[s*np*dim], np*dim);CHKERRQ(ierr); 1835c0ce001eSMatthew G. Knepley ierr = PetscArraycpy(w, &weights[s*np*Nc], np*Nc);CHKERRQ(ierr); 1836c0ce001eSMatthew G. Knepley ierr = PetscQuadratureSetData(qs, dim, Nc, np, p, w);CHKERRQ(ierr); 1837c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceSimpleSetFunctional(Qref, s, qs);CHKERRQ(ierr); 1838c0ce001eSMatthew G. Knepley ierr = PetscQuadratureDestroy(&qs);CHKERRQ(ierr); 1839c0ce001eSMatthew G. Knepley } 1840c0ce001eSMatthew G. Knepley ierr = PetscFVSetQuadrature(*fvRef, qref);CHKERRQ(ierr); 1841*412e9a14SMatthew G. Knepley ierr = DMPlexCellRefinerDestroy(&cr);CHKERRQ(ierr); 1842c0ce001eSMatthew G. Knepley ierr = PetscQuadratureDestroy(&qref);CHKERRQ(ierr); 1843c0ce001eSMatthew G. Knepley ierr = PetscDualSpaceDestroy(&Qref);CHKERRQ(ierr); 1844c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1845c0ce001eSMatthew G. Knepley } 1846c0ce001eSMatthew G. Knepley 184788f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVDestroy_Upwind(PetscFV fvm) 1848c0ce001eSMatthew G. Knepley { 1849c0ce001eSMatthew G. Knepley PetscFV_Upwind *b = (PetscFV_Upwind *) fvm->data; 1850c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1851c0ce001eSMatthew G. Knepley 1852c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1853c0ce001eSMatthew G. Knepley ierr = PetscFree(b);CHKERRQ(ierr); 1854c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1855c0ce001eSMatthew G. Knepley } 1856c0ce001eSMatthew G. Knepley 185788f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVView_Upwind_Ascii(PetscFV fv, PetscViewer viewer) 1858c0ce001eSMatthew G. Knepley { 1859c0ce001eSMatthew G. Knepley PetscInt Nc, c; 1860c0ce001eSMatthew G. Knepley PetscViewerFormat format; 1861c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1862c0ce001eSMatthew G. Knepley 1863c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1864c0ce001eSMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 1865c0ce001eSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 1866c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Upwind Finite Volume:\n");CHKERRQ(ierr); 1867c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, " num components: %d\n", Nc);CHKERRQ(ierr); 1868c0ce001eSMatthew G. Knepley for (c = 0; c < Nc; c++) { 1869c0ce001eSMatthew G. Knepley if (fv->componentNames[c]) { 1870c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, " component %d: %s\n", c, fv->componentNames[c]);CHKERRQ(ierr); 1871c0ce001eSMatthew G. Knepley } 1872c0ce001eSMatthew G. Knepley } 1873c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1874c0ce001eSMatthew G. Knepley } 1875c0ce001eSMatthew G. Knepley 187688f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVView_Upwind(PetscFV fv, PetscViewer viewer) 1877c0ce001eSMatthew G. Knepley { 1878c0ce001eSMatthew G. Knepley PetscBool iascii; 1879c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1880c0ce001eSMatthew G. Knepley 1881c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1882c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fv, PETSCFV_CLASSID, 1); 1883c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1884c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 1885c0ce001eSMatthew G. Knepley if (iascii) {ierr = PetscFVView_Upwind_Ascii(fv, viewer);CHKERRQ(ierr);} 1886c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1887c0ce001eSMatthew G. Knepley } 1888c0ce001eSMatthew G. Knepley 188988f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVSetUp_Upwind(PetscFV fvm) 1890c0ce001eSMatthew G. Knepley { 1891c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1892c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1893c0ce001eSMatthew G. Knepley } 1894c0ce001eSMatthew G. Knepley 1895c0ce001eSMatthew G. Knepley /* 1896c0ce001eSMatthew G. Knepley neighborVol[f*2+0] contains the left geom 1897c0ce001eSMatthew G. Knepley neighborVol[f*2+1] contains the right geom 1898c0ce001eSMatthew G. Knepley */ 189988f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVIntegrateRHSFunction_Upwind(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol, 1900c0ce001eSMatthew G. Knepley PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[]) 1901c0ce001eSMatthew G. Knepley { 1902c0ce001eSMatthew G. Knepley void (*riemann)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *); 1903c0ce001eSMatthew G. Knepley void *rctx; 1904c0ce001eSMatthew G. Knepley PetscScalar *flux = fvm->fluxWork; 1905c0ce001eSMatthew G. Knepley const PetscScalar *constants; 1906c0ce001eSMatthew G. Knepley PetscInt dim, numConstants, pdim, totDim, Nc, off, f, d; 1907c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1908c0ce001eSMatthew G. Knepley 1909c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1910c0ce001eSMatthew G. Knepley ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr); 1911c0ce001eSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1912c0ce001eSMatthew G. Knepley ierr = PetscDSGetFieldOffset(prob, field, &off);CHKERRQ(ierr); 1913c0ce001eSMatthew G. Knepley ierr = PetscDSGetRiemannSolver(prob, field, &riemann);CHKERRQ(ierr); 1914c0ce001eSMatthew G. Knepley ierr = PetscDSGetContext(prob, field, &rctx);CHKERRQ(ierr); 1915c0ce001eSMatthew G. Knepley ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 1916c0ce001eSMatthew G. Knepley ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr); 1917c0ce001eSMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 1918c0ce001eSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 1919c0ce001eSMatthew G. Knepley (*riemann)(dim, pdim, fgeom[f].centroid, fgeom[f].normal, &uL[f*Nc], &uR[f*Nc], numConstants, constants, flux, rctx); 1920c0ce001eSMatthew G. Knepley for (d = 0; d < pdim; ++d) { 1921c0ce001eSMatthew G. Knepley fluxL[f*totDim+off+d] = flux[d] / neighborVol[f*2+0]; 1922c0ce001eSMatthew G. Knepley fluxR[f*totDim+off+d] = flux[d] / neighborVol[f*2+1]; 1923c0ce001eSMatthew G. Knepley } 1924c0ce001eSMatthew G. Knepley } 1925c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1926c0ce001eSMatthew G. Knepley } 1927c0ce001eSMatthew G. Knepley 192888f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVInitialize_Upwind(PetscFV fvm) 1929c0ce001eSMatthew G. Knepley { 1930c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1931c0ce001eSMatthew G. Knepley fvm->ops->setfromoptions = NULL; 1932c0ce001eSMatthew G. Knepley fvm->ops->setup = PetscFVSetUp_Upwind; 1933c0ce001eSMatthew G. Knepley fvm->ops->view = PetscFVView_Upwind; 1934c0ce001eSMatthew G. Knepley fvm->ops->destroy = PetscFVDestroy_Upwind; 1935c0ce001eSMatthew G. Knepley fvm->ops->integraterhsfunction = PetscFVIntegrateRHSFunction_Upwind; 1936c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1937c0ce001eSMatthew G. Knepley } 1938c0ce001eSMatthew G. Knepley 1939c0ce001eSMatthew G. Knepley /*MC 1940c0ce001eSMatthew G. Knepley PETSCFVUPWIND = "upwind" - A PetscFV object 1941c0ce001eSMatthew G. Knepley 1942c0ce001eSMatthew G. Knepley Level: intermediate 1943c0ce001eSMatthew G. Knepley 1944c0ce001eSMatthew G. Knepley .seealso: PetscFVType, PetscFVCreate(), PetscFVSetType() 1945c0ce001eSMatthew G. Knepley M*/ 1946c0ce001eSMatthew G. Knepley 1947c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFVCreate_Upwind(PetscFV fvm) 1948c0ce001eSMatthew G. Knepley { 1949c0ce001eSMatthew G. Knepley PetscFV_Upwind *b; 1950c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1951c0ce001eSMatthew G. Knepley 1952c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1953c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1954c0ce001eSMatthew G. Knepley ierr = PetscNewLog(fvm,&b);CHKERRQ(ierr); 1955c0ce001eSMatthew G. Knepley fvm->data = b; 1956c0ce001eSMatthew G. Knepley 1957c0ce001eSMatthew G. Knepley ierr = PetscFVInitialize_Upwind(fvm);CHKERRQ(ierr); 1958c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1959c0ce001eSMatthew G. Knepley } 1960c0ce001eSMatthew G. Knepley 1961c0ce001eSMatthew G. Knepley #include <petscblaslapack.h> 1962c0ce001eSMatthew G. Knepley 196388f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVDestroy_LeastSquares(PetscFV fvm) 1964c0ce001eSMatthew G. Knepley { 1965c0ce001eSMatthew G. Knepley PetscFV_LeastSquares *ls = (PetscFV_LeastSquares *) fvm->data; 1966c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1967c0ce001eSMatthew G. Knepley 1968c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1969c0ce001eSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject) fvm, "PetscFVLeastSquaresSetMaxFaces_C", NULL);CHKERRQ(ierr); 1970c0ce001eSMatthew G. Knepley ierr = PetscFree4(ls->B, ls->Binv, ls->tau, ls->work);CHKERRQ(ierr); 1971c0ce001eSMatthew G. Knepley ierr = PetscFree(ls);CHKERRQ(ierr); 1972c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1973c0ce001eSMatthew G. Knepley } 1974c0ce001eSMatthew G. Knepley 197588f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVView_LeastSquares_Ascii(PetscFV fv, PetscViewer viewer) 1976c0ce001eSMatthew G. Knepley { 1977c0ce001eSMatthew G. Knepley PetscInt Nc, c; 1978c0ce001eSMatthew G. Knepley PetscViewerFormat format; 1979c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1980c0ce001eSMatthew G. Knepley 1981c0ce001eSMatthew G. Knepley PetscFunctionBegin; 1982c0ce001eSMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 1983c0ce001eSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 1984c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Finite Volume with Least Squares Reconstruction:\n");CHKERRQ(ierr); 1985c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, " num components: %d\n", Nc);CHKERRQ(ierr); 1986c0ce001eSMatthew G. Knepley for (c = 0; c < Nc; c++) { 1987c0ce001eSMatthew G. Knepley if (fv->componentNames[c]) { 1988c0ce001eSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, " component %d: %s\n", c, fv->componentNames[c]);CHKERRQ(ierr); 1989c0ce001eSMatthew G. Knepley } 1990c0ce001eSMatthew G. Knepley } 1991c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 1992c0ce001eSMatthew G. Knepley } 1993c0ce001eSMatthew G. Knepley 199488f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVView_LeastSquares(PetscFV fv, PetscViewer viewer) 1995c0ce001eSMatthew G. Knepley { 1996c0ce001eSMatthew G. Knepley PetscBool iascii; 1997c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 1998c0ce001eSMatthew G. Knepley 1999c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2000c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fv, PETSCFV_CLASSID, 1); 2001c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2002c0ce001eSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 2003c0ce001eSMatthew G. Knepley if (iascii) {ierr = PetscFVView_LeastSquares_Ascii(fv, viewer);CHKERRQ(ierr);} 2004c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2005c0ce001eSMatthew G. Knepley } 2006c0ce001eSMatthew G. Knepley 200788f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVSetUp_LeastSquares(PetscFV fvm) 2008c0ce001eSMatthew G. Knepley { 2009c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2010c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2011c0ce001eSMatthew G. Knepley } 2012c0ce001eSMatthew G. Knepley 2013c0ce001eSMatthew G. Knepley /* Overwrites A. Can only handle full-rank problems with m>=n */ 2014c0ce001eSMatthew G. Knepley static PetscErrorCode PetscFVLeastSquaresPseudoInverse_Static(PetscInt m,PetscInt mstride,PetscInt n,PetscScalar *A,PetscScalar *Ainv,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 2015c0ce001eSMatthew G. Knepley { 2016c0ce001eSMatthew G. Knepley PetscBool debug = PETSC_FALSE; 2017c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 2018c0ce001eSMatthew G. Knepley PetscBLASInt M,N,K,lda,ldb,ldwork,info; 2019c0ce001eSMatthew G. Knepley PetscScalar *R,*Q,*Aback,Alpha; 2020c0ce001eSMatthew G. Knepley 2021c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2022c0ce001eSMatthew G. Knepley if (debug) { 2023c0ce001eSMatthew G. Knepley ierr = PetscMalloc1(m*n,&Aback);CHKERRQ(ierr); 2024c0ce001eSMatthew G. Knepley ierr = PetscArraycpy(Aback,A,m*n);CHKERRQ(ierr); 2025c0ce001eSMatthew G. Knepley } 2026c0ce001eSMatthew G. Knepley 2027c0ce001eSMatthew G. Knepley ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 2028c0ce001eSMatthew G. Knepley ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 2029c0ce001eSMatthew G. Knepley ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 2030c0ce001eSMatthew G. Knepley ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 2031c0ce001eSMatthew G. Knepley ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 2032c0ce001eSMatthew G. Knepley LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info); 2033c0ce001eSMatthew G. Knepley ierr = PetscFPTrapPop();CHKERRQ(ierr); 2034c0ce001eSMatthew G. Knepley if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 2035c0ce001eSMatthew G. Knepley R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 2036c0ce001eSMatthew G. Knepley 2037c0ce001eSMatthew G. Knepley /* Extract an explicit representation of Q */ 2038c0ce001eSMatthew G. Knepley Q = Ainv; 2039c0ce001eSMatthew G. Knepley ierr = PetscArraycpy(Q,A,mstride*n);CHKERRQ(ierr); 2040c0ce001eSMatthew G. Knepley K = N; /* full rank */ 2041c0ce001eSMatthew G. Knepley PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 2042c0ce001eSMatthew G. Knepley if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 2043c0ce001eSMatthew G. Knepley 2044c0ce001eSMatthew G. Knepley /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 2045c0ce001eSMatthew G. Knepley Alpha = 1.0; 2046c0ce001eSMatthew G. Knepley ldb = lda; 2047c0ce001eSMatthew G. Knepley BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb); 2048c0ce001eSMatthew G. Knepley /* Ainv is Q, overwritten with inverse */ 2049c0ce001eSMatthew G. Knepley 2050c0ce001eSMatthew G. Knepley if (debug) { /* Check that pseudo-inverse worked */ 2051c0ce001eSMatthew G. Knepley PetscScalar Beta = 0.0; 2052c0ce001eSMatthew G. Knepley PetscBLASInt ldc; 2053c0ce001eSMatthew G. Knepley K = N; 2054c0ce001eSMatthew G. Knepley ldc = N; 2055c0ce001eSMatthew G. Knepley BLASgemm_("ConjugateTranspose","Normal",&N,&K,&M,&Alpha,Ainv,&lda,Aback,&ldb,&Beta,work,&ldc); 2056c0ce001eSMatthew G. Knepley ierr = PetscScalarView(n*n,work,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 2057c0ce001eSMatthew G. Knepley ierr = PetscFree(Aback);CHKERRQ(ierr); 2058c0ce001eSMatthew G. Knepley } 2059c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2060c0ce001eSMatthew G. Knepley } 2061c0ce001eSMatthew G. Knepley 2062c0ce001eSMatthew G. Knepley /* Overwrites A. Can handle degenerate problems and m<n. */ 2063c0ce001eSMatthew G. Knepley static PetscErrorCode PetscFVLeastSquaresPseudoInverseSVD_Static(PetscInt m,PetscInt mstride,PetscInt n,PetscScalar *A,PetscScalar *Ainv,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 2064c0ce001eSMatthew G. Knepley { 2065c0ce001eSMatthew G. Knepley PetscBool debug = PETSC_FALSE; 2066c0ce001eSMatthew G. Knepley PetscScalar *Brhs, *Aback; 2067c0ce001eSMatthew G. Knepley PetscScalar *tmpwork; 2068c0ce001eSMatthew G. Knepley PetscReal rcond; 2069c0ce001eSMatthew G. Knepley #if defined (PETSC_USE_COMPLEX) 2070c0ce001eSMatthew G. Knepley PetscInt rworkSize; 2071c0ce001eSMatthew G. Knepley PetscReal *rwork; 2072c0ce001eSMatthew G. Knepley #endif 2073c0ce001eSMatthew G. Knepley PetscInt i, j, maxmn; 2074c0ce001eSMatthew G. Knepley PetscBLASInt M, N, lda, ldb, ldwork; 2075c0ce001eSMatthew G. Knepley PetscBLASInt nrhs, irank, info; 2076c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 2077c0ce001eSMatthew G. Knepley 2078c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2079c0ce001eSMatthew G. Knepley if (debug) { 2080c0ce001eSMatthew G. Knepley ierr = PetscMalloc1(m*n,&Aback);CHKERRQ(ierr); 2081c0ce001eSMatthew G. Knepley ierr = PetscArraycpy(Aback,A,m*n);CHKERRQ(ierr); 2082c0ce001eSMatthew G. Knepley } 2083c0ce001eSMatthew G. Knepley 2084c0ce001eSMatthew G. Knepley /* initialize to identity */ 2085c0ce001eSMatthew G. Knepley tmpwork = Ainv; 2086c0ce001eSMatthew G. Knepley Brhs = work; 2087c0ce001eSMatthew G. Knepley maxmn = PetscMax(m,n); 2088c0ce001eSMatthew G. Knepley for (j=0; j<maxmn; j++) { 2089c0ce001eSMatthew G. Knepley for (i=0; i<maxmn; i++) Brhs[i + j*maxmn] = 1.0*(i == j); 2090c0ce001eSMatthew G. Knepley } 2091c0ce001eSMatthew G. Knepley 2092c0ce001eSMatthew G. Knepley ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 2093c0ce001eSMatthew G. Knepley ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 2094c0ce001eSMatthew G. Knepley ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 2095c0ce001eSMatthew G. Knepley ierr = PetscBLASIntCast(maxmn,&ldb);CHKERRQ(ierr); 2096c0ce001eSMatthew G. Knepley ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 2097c0ce001eSMatthew G. Knepley rcond = -1; 2098c0ce001eSMatthew G. Knepley ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 2099c0ce001eSMatthew G. Knepley nrhs = M; 2100c0ce001eSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX) 2101c0ce001eSMatthew G. Knepley rworkSize = 5 * PetscMin(M,N); 2102c0ce001eSMatthew G. Knepley ierr = PetscMalloc1(rworkSize,&rwork);CHKERRQ(ierr); 2103c0ce001eSMatthew G. Knepley LAPACKgelss_(&M,&N,&nrhs,A,&lda,Brhs,&ldb, (PetscReal *) tau,&rcond,&irank,tmpwork,&ldwork,rwork,&info); 2104c0ce001eSMatthew G. Knepley ierr = PetscFPTrapPop();CHKERRQ(ierr); 2105c0ce001eSMatthew G. Knepley ierr = PetscFree(rwork);CHKERRQ(ierr); 2106c0ce001eSMatthew G. Knepley #else 2107c0ce001eSMatthew G. Knepley nrhs = M; 2108c0ce001eSMatthew G. Knepley LAPACKgelss_(&M,&N,&nrhs,A,&lda,Brhs,&ldb, (PetscReal *) tau,&rcond,&irank,tmpwork,&ldwork,&info); 2109c0ce001eSMatthew G. Knepley ierr = PetscFPTrapPop();CHKERRQ(ierr); 2110c0ce001eSMatthew G. Knepley #endif 2111c0ce001eSMatthew G. Knepley if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGELSS error"); 2112c0ce001eSMatthew G. Knepley /* The following check should be turned into a diagnostic as soon as someone wants to do this intentionally */ 2113c0ce001eSMatthew G. Knepley if (irank < PetscMin(M,N)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Rank deficient least squares fit, indicates an isolated cell with two colinear points"); 2114c0ce001eSMatthew G. Knepley 2115c0ce001eSMatthew G. Knepley /* Brhs shaped (M,nrhs) column-major coldim=mstride was overwritten by Ainv shaped (N,nrhs) column-major coldim=maxmn. 2116c0ce001eSMatthew G. Knepley * Here we transpose to (N,nrhs) row-major rowdim=mstride. */ 2117c0ce001eSMatthew G. Knepley for (i=0; i<n; i++) { 2118c0ce001eSMatthew G. Knepley for (j=0; j<nrhs; j++) Ainv[i*mstride+j] = Brhs[i + j*maxmn]; 2119c0ce001eSMatthew G. Knepley } 2120c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2121c0ce001eSMatthew G. Knepley } 2122c0ce001eSMatthew G. Knepley 2123c0ce001eSMatthew G. Knepley #if 0 2124c0ce001eSMatthew G. Knepley static PetscErrorCode PetscFVLeastSquaresDebugCell_Static(PetscFV fvm, PetscInt cell, DM dm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 2125c0ce001eSMatthew G. Knepley { 2126c0ce001eSMatthew G. Knepley PetscReal grad[2] = {0, 0}; 2127c0ce001eSMatthew G. Knepley const PetscInt *faces; 2128c0ce001eSMatthew G. Knepley PetscInt numFaces, f; 2129c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 2130c0ce001eSMatthew G. Knepley 2131c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2132c0ce001eSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 2133c0ce001eSMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 2134c0ce001eSMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 2135c0ce001eSMatthew G. Knepley const PetscInt *fcells; 2136c0ce001eSMatthew G. Knepley const CellGeom *cg1; 2137c0ce001eSMatthew G. Knepley const FaceGeom *fg; 2138c0ce001eSMatthew G. Knepley 2139c0ce001eSMatthew G. Knepley ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 2140c0ce001eSMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 2141c0ce001eSMatthew G. Knepley for (i = 0; i < 2; ++i) { 2142c0ce001eSMatthew G. Knepley PetscScalar du; 2143c0ce001eSMatthew G. Knepley 2144c0ce001eSMatthew G. Knepley if (fcells[i] == c) continue; 2145c0ce001eSMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, fcells[i], cgeom, &cg1);CHKERRQ(ierr); 2146c0ce001eSMatthew G. Knepley du = cg1->centroid[0] + 3*cg1->centroid[1] - (cg->centroid[0] + 3*cg->centroid[1]); 2147c0ce001eSMatthew G. Knepley grad[0] += fg->grad[!i][0] * du; 2148c0ce001eSMatthew G. Knepley grad[1] += fg->grad[!i][1] * du; 2149c0ce001eSMatthew G. Knepley } 2150c0ce001eSMatthew G. Knepley } 2151c0ce001eSMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "cell[%d] grad (%g, %g)\n", cell, grad[0], grad[1]);CHKERRQ(ierr); 2152c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2153c0ce001eSMatthew G. Knepley } 2154c0ce001eSMatthew G. Knepley #endif 2155c0ce001eSMatthew G. Knepley 2156c0ce001eSMatthew G. Knepley /* 2157c0ce001eSMatthew G. Knepley PetscFVComputeGradient - Compute the gradient reconstruction matrix for a given cell 2158c0ce001eSMatthew G. Knepley 2159c0ce001eSMatthew G. Knepley Input Parameters: 2160c0ce001eSMatthew G. Knepley + fvm - The PetscFV object 2161c0ce001eSMatthew G. Knepley . numFaces - The number of cell faces which are not constrained 2162c0ce001eSMatthew G. Knepley . dx - The vector from the cell centroid to the neighboring cell centroid for each face 2163c0ce001eSMatthew G. Knepley 2164c0ce001eSMatthew G. Knepley Level: developer 2165c0ce001eSMatthew G. Knepley 2166c0ce001eSMatthew G. Knepley .seealso: PetscFVCreate() 2167c0ce001eSMatthew G. Knepley */ 216888f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVComputeGradient_LeastSquares(PetscFV fvm, PetscInt numFaces, const PetscScalar dx[], PetscScalar grad[]) 2169c0ce001eSMatthew G. Knepley { 2170c0ce001eSMatthew G. Knepley PetscFV_LeastSquares *ls = (PetscFV_LeastSquares *) fvm->data; 2171c0ce001eSMatthew G. Knepley const PetscBool useSVD = PETSC_TRUE; 2172c0ce001eSMatthew G. Knepley const PetscInt maxFaces = ls->maxFaces; 2173c0ce001eSMatthew G. Knepley PetscInt dim, f, d; 2174c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 2175c0ce001eSMatthew G. Knepley 2176c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2177c0ce001eSMatthew G. Knepley if (numFaces > maxFaces) { 2178c0ce001eSMatthew G. Knepley if (maxFaces < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Reconstruction has not been initialized, call PetscFVLeastSquaresSetMaxFaces()"); 2179c0ce001eSMatthew G. Knepley SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input faces %D > %D maxfaces", numFaces, maxFaces); 2180c0ce001eSMatthew G. Knepley } 2181c0ce001eSMatthew G. Knepley ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr); 2182c0ce001eSMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 2183c0ce001eSMatthew G. Knepley for (d = 0; d < dim; ++d) ls->B[d*maxFaces+f] = dx[f*dim+d]; 2184c0ce001eSMatthew G. Knepley } 2185c0ce001eSMatthew G. Knepley /* Overwrites B with garbage, returns Binv in row-major format */ 2186c0ce001eSMatthew G. Knepley if (useSVD) {ierr = PetscFVLeastSquaresPseudoInverseSVD_Static(numFaces, maxFaces, dim, ls->B, ls->Binv, ls->tau, ls->workSize, ls->work);CHKERRQ(ierr);} 2187c0ce001eSMatthew G. Knepley else {ierr = PetscFVLeastSquaresPseudoInverse_Static(numFaces, maxFaces, dim, ls->B, ls->Binv, ls->tau, ls->workSize, ls->work);CHKERRQ(ierr);} 2188c0ce001eSMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 2189c0ce001eSMatthew G. Knepley for (d = 0; d < dim; ++d) grad[f*dim+d] = ls->Binv[d*maxFaces+f]; 2190c0ce001eSMatthew G. Knepley } 2191c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2192c0ce001eSMatthew G. Knepley } 2193c0ce001eSMatthew G. Knepley 2194c0ce001eSMatthew G. Knepley /* 2195c0ce001eSMatthew G. Knepley neighborVol[f*2+0] contains the left geom 2196c0ce001eSMatthew G. Knepley neighborVol[f*2+1] contains the right geom 2197c0ce001eSMatthew G. Knepley */ 219888f5f89eSMatthew G. Knepley static PetscErrorCode PetscFVIntegrateRHSFunction_LeastSquares(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol, 2199c0ce001eSMatthew G. Knepley PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[]) 2200c0ce001eSMatthew G. Knepley { 2201c0ce001eSMatthew G. Knepley void (*riemann)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *); 2202c0ce001eSMatthew G. Knepley void *rctx; 2203c0ce001eSMatthew G. Knepley PetscScalar *flux = fvm->fluxWork; 2204c0ce001eSMatthew G. Knepley const PetscScalar *constants; 2205c0ce001eSMatthew G. Knepley PetscInt dim, numConstants, pdim, Nc, totDim, off, f, d; 2206c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 2207c0ce001eSMatthew G. Knepley 2208c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2209c0ce001eSMatthew G. Knepley ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr); 2210c0ce001eSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 2211c0ce001eSMatthew G. Knepley ierr = PetscDSGetFieldOffset(prob, field, &off);CHKERRQ(ierr); 2212c0ce001eSMatthew G. Knepley ierr = PetscDSGetRiemannSolver(prob, field, &riemann);CHKERRQ(ierr); 2213c0ce001eSMatthew G. Knepley ierr = PetscDSGetContext(prob, field, &rctx);CHKERRQ(ierr); 2214c0ce001eSMatthew G. Knepley ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 2215c0ce001eSMatthew G. Knepley ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr); 2216c0ce001eSMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 2217c0ce001eSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 2218c0ce001eSMatthew G. Knepley (*riemann)(dim, pdim, fgeom[f].centroid, fgeom[f].normal, &uL[f*Nc], &uR[f*Nc], numConstants, constants, flux, rctx); 2219c0ce001eSMatthew G. Knepley for (d = 0; d < pdim; ++d) { 2220c0ce001eSMatthew G. Knepley fluxL[f*totDim+off+d] = flux[d] / neighborVol[f*2+0]; 2221c0ce001eSMatthew G. Knepley fluxR[f*totDim+off+d] = flux[d] / neighborVol[f*2+1]; 2222c0ce001eSMatthew G. Knepley } 2223c0ce001eSMatthew G. Knepley } 2224c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2225c0ce001eSMatthew G. Knepley } 2226c0ce001eSMatthew G. Knepley 2227c0ce001eSMatthew G. Knepley static PetscErrorCode PetscFVLeastSquaresSetMaxFaces_LS(PetscFV fvm, PetscInt maxFaces) 2228c0ce001eSMatthew G. Knepley { 2229c0ce001eSMatthew G. Knepley PetscFV_LeastSquares *ls = (PetscFV_LeastSquares *) fvm->data; 2230c0ce001eSMatthew G. Knepley PetscInt dim, m, n, nrhs, minwork; 2231c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 2232c0ce001eSMatthew G. Knepley 2233c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2234c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 2235c0ce001eSMatthew G. Knepley ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr); 2236c0ce001eSMatthew G. Knepley ierr = PetscFree4(ls->B, ls->Binv, ls->tau, ls->work);CHKERRQ(ierr); 2237c0ce001eSMatthew G. Knepley ls->maxFaces = maxFaces; 2238c0ce001eSMatthew G. Knepley m = ls->maxFaces; 2239c0ce001eSMatthew G. Knepley n = dim; 2240c0ce001eSMatthew G. Knepley nrhs = ls->maxFaces; 2241c0ce001eSMatthew G. Knepley minwork = 3*PetscMin(m,n) + PetscMax(2*PetscMin(m,n), PetscMax(PetscMax(m,n), nrhs)); /* required by LAPACK */ 2242c0ce001eSMatthew G. Knepley ls->workSize = 5*minwork; /* We can afford to be extra generous */ 2243c0ce001eSMatthew G. Knepley ierr = PetscMalloc4(ls->maxFaces*dim,&ls->B,ls->workSize,&ls->Binv,ls->maxFaces,&ls->tau,ls->workSize,&ls->work);CHKERRQ(ierr); 2244c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2245c0ce001eSMatthew G. Knepley } 2246c0ce001eSMatthew G. Knepley 2247c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVInitialize_LeastSquares(PetscFV fvm) 2248c0ce001eSMatthew G. Knepley { 2249c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2250c0ce001eSMatthew G. Knepley fvm->ops->setfromoptions = NULL; 2251c0ce001eSMatthew G. Knepley fvm->ops->setup = PetscFVSetUp_LeastSquares; 2252c0ce001eSMatthew G. Knepley fvm->ops->view = PetscFVView_LeastSquares; 2253c0ce001eSMatthew G. Knepley fvm->ops->destroy = PetscFVDestroy_LeastSquares; 2254c0ce001eSMatthew G. Knepley fvm->ops->computegradient = PetscFVComputeGradient_LeastSquares; 2255c0ce001eSMatthew G. Knepley fvm->ops->integraterhsfunction = PetscFVIntegrateRHSFunction_LeastSquares; 2256c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2257c0ce001eSMatthew G. Knepley } 2258c0ce001eSMatthew G. Knepley 2259c0ce001eSMatthew G. Knepley /*MC 2260c0ce001eSMatthew G. Knepley PETSCFVLEASTSQUARES = "leastsquares" - A PetscFV object 2261c0ce001eSMatthew G. Knepley 2262c0ce001eSMatthew G. Knepley Level: intermediate 2263c0ce001eSMatthew G. Knepley 2264c0ce001eSMatthew G. Knepley .seealso: PetscFVType, PetscFVCreate(), PetscFVSetType() 2265c0ce001eSMatthew G. Knepley M*/ 2266c0ce001eSMatthew G. Knepley 2267c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFVCreate_LeastSquares(PetscFV fvm) 2268c0ce001eSMatthew G. Knepley { 2269c0ce001eSMatthew G. Knepley PetscFV_LeastSquares *ls; 2270c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 2271c0ce001eSMatthew G. Knepley 2272c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2273c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 2274c0ce001eSMatthew G. Knepley ierr = PetscNewLog(fvm, &ls);CHKERRQ(ierr); 2275c0ce001eSMatthew G. Knepley fvm->data = ls; 2276c0ce001eSMatthew G. Knepley 2277c0ce001eSMatthew G. Knepley ls->maxFaces = -1; 2278c0ce001eSMatthew G. Knepley ls->workSize = -1; 2279c0ce001eSMatthew G. Knepley ls->B = NULL; 2280c0ce001eSMatthew G. Knepley ls->Binv = NULL; 2281c0ce001eSMatthew G. Knepley ls->tau = NULL; 2282c0ce001eSMatthew G. Knepley ls->work = NULL; 2283c0ce001eSMatthew G. Knepley 2284c0ce001eSMatthew G. Knepley ierr = PetscFVSetComputeGradients(fvm, PETSC_TRUE);CHKERRQ(ierr); 2285c0ce001eSMatthew G. Knepley ierr = PetscFVInitialize_LeastSquares(fvm);CHKERRQ(ierr); 2286c0ce001eSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject) fvm, "PetscFVLeastSquaresSetMaxFaces_C", PetscFVLeastSquaresSetMaxFaces_LS);CHKERRQ(ierr); 2287c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2288c0ce001eSMatthew G. Knepley } 2289c0ce001eSMatthew G. Knepley 2290c0ce001eSMatthew G. Knepley /*@ 2291c0ce001eSMatthew G. Knepley PetscFVLeastSquaresSetMaxFaces - Set the maximum number of cell faces for gradient reconstruction 2292c0ce001eSMatthew G. Knepley 2293c0ce001eSMatthew G. Knepley Not collective 2294c0ce001eSMatthew G. Knepley 2295c0ce001eSMatthew G. Knepley Input parameters: 2296c0ce001eSMatthew G. Knepley + fvm - The PetscFV object 2297c0ce001eSMatthew G. Knepley - maxFaces - The maximum number of cell faces 2298c0ce001eSMatthew G. Knepley 2299c0ce001eSMatthew G. Knepley Level: intermediate 2300c0ce001eSMatthew G. Knepley 2301c0ce001eSMatthew G. Knepley .seealso: PetscFVCreate(), PETSCFVLEASTSQUARES 2302c0ce001eSMatthew G. Knepley @*/ 2303c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVLeastSquaresSetMaxFaces(PetscFV fvm, PetscInt maxFaces) 2304c0ce001eSMatthew G. Knepley { 2305c0ce001eSMatthew G. Knepley PetscErrorCode ierr; 2306c0ce001eSMatthew G. Knepley 2307c0ce001eSMatthew G. Knepley PetscFunctionBegin; 2308c0ce001eSMatthew G. Knepley PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 2309c0ce001eSMatthew G. Knepley ierr = PetscTryMethod(fvm, "PetscFVLeastSquaresSetMaxFaces_C", (PetscFV,PetscInt), (fvm,maxFaces));CHKERRQ(ierr); 2310c0ce001eSMatthew G. Knepley PetscFunctionReturn(0); 2311c0ce001eSMatthew G. Knepley } 2312