1 #include <petsc-private/petscdsimpl.h> /*I "petscds.h" I*/ 2 3 PetscClassId PETSCDS_CLASSID = 0; 4 5 PetscFunctionList PetscDSList = NULL; 6 PetscBool PetscDSRegisterAllCalled = PETSC_FALSE; 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "PetscDSRegister" 10 /*@C 11 PetscDSRegister - Adds a new PetscDS implementation 12 13 Not Collective 14 15 Input Parameters: 16 + name - The name of a new user-defined creation routine 17 - create_func - The creation routine itself 18 19 Notes: 20 PetscDSRegister() may be called multiple times to add several user-defined PetscDSs 21 22 Sample usage: 23 .vb 24 PetscDSRegister("my_ds", MyPetscDSCreate); 25 .ve 26 27 Then, your PetscDS type can be chosen with the procedural interface via 28 .vb 29 PetscDSCreate(MPI_Comm, PetscDS *); 30 PetscDSSetType(PetscDS, "my_ds"); 31 .ve 32 or at runtime via the option 33 .vb 34 -petscds_type my_ds 35 .ve 36 37 Level: advanced 38 39 .keywords: PetscDS, register 40 .seealso: PetscDSRegisterAll(), PetscDSRegisterDestroy() 41 42 @*/ 43 PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS)) 44 { 45 PetscErrorCode ierr; 46 47 PetscFunctionBegin; 48 ierr = PetscFunctionListAdd(&PetscDSList, sname, function);CHKERRQ(ierr); 49 PetscFunctionReturn(0); 50 } 51 52 #undef __FUNCT__ 53 #define __FUNCT__ "PetscDSSetType" 54 /*@C 55 PetscDSSetType - Builds a particular PetscDS 56 57 Collective on PetscDS 58 59 Input Parameters: 60 + prob - The PetscDS object 61 - name - The kind of system 62 63 Options Database Key: 64 . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types 65 66 Level: intermediate 67 68 .keywords: PetscDS, set, type 69 .seealso: PetscDSGetType(), PetscDSCreate() 70 @*/ 71 PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name) 72 { 73 PetscErrorCode (*r)(PetscDS); 74 PetscBool match; 75 PetscErrorCode ierr; 76 77 PetscFunctionBegin; 78 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 79 ierr = PetscObjectTypeCompare((PetscObject) prob, name, &match);CHKERRQ(ierr); 80 if (match) PetscFunctionReturn(0); 81 82 ierr = PetscDSRegisterAll();CHKERRQ(ierr); 83 ierr = PetscFunctionListFind(PetscDSList, name, &r);CHKERRQ(ierr); 84 if (!r) SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name); 85 86 if (prob->ops->destroy) { 87 ierr = (*prob->ops->destroy)(prob);CHKERRQ(ierr); 88 prob->ops->destroy = NULL; 89 } 90 ierr = (*r)(prob);CHKERRQ(ierr); 91 ierr = PetscObjectChangeTypeName((PetscObject) prob, name);CHKERRQ(ierr); 92 PetscFunctionReturn(0); 93 } 94 95 #undef __FUNCT__ 96 #define __FUNCT__ "PetscDSGetType" 97 /*@C 98 PetscDSGetType - Gets the PetscDS type name (as a string) from the object. 99 100 Not Collective 101 102 Input Parameter: 103 . prob - The PetscDS 104 105 Output Parameter: 106 . name - The PetscDS type name 107 108 Level: intermediate 109 110 .keywords: PetscDS, get, type, name 111 .seealso: PetscDSSetType(), PetscDSCreate() 112 @*/ 113 PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name) 114 { 115 PetscErrorCode ierr; 116 117 PetscFunctionBegin; 118 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 119 PetscValidPointer(name, 2); 120 ierr = PetscDSRegisterAll();CHKERRQ(ierr); 121 *name = ((PetscObject) prob)->type_name; 122 PetscFunctionReturn(0); 123 } 124 125 #undef __FUNCT__ 126 #define __FUNCT__ "PetscDSView_Ascii" 127 static PetscErrorCode PetscDSView_Ascii(PetscDS prob, PetscViewer viewer) 128 { 129 PetscViewerFormat format; 130 PetscInt f; 131 PetscErrorCode ierr; 132 133 PetscFunctionBegin; 134 ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 135 ierr = PetscViewerASCIIPrintf(viewer, "Discrete System with %d fields\n", prob->Nf);CHKERRQ(ierr); 136 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 137 for (f = 0; f < prob->Nf; ++f) { 138 PetscObject obj; 139 PetscClassId id; 140 const char *name; 141 PetscInt Nc; 142 143 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 144 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 145 ierr = PetscObjectGetName(obj, &name);CHKERRQ(ierr); 146 ierr = PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>");CHKERRQ(ierr); 147 if (id == PETSCFE_CLASSID) { 148 ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr); 149 ierr = PetscViewerASCIIPrintf(viewer, " FEM");CHKERRQ(ierr); 150 } else if (id == PETSCFV_CLASSID) { 151 ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr); 152 ierr = PetscViewerASCIIPrintf(viewer, " FVM");CHKERRQ(ierr); 153 } 154 else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 155 if (Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "%d components", Nc);CHKERRQ(ierr);} 156 else {ierr = PetscViewerASCIIPrintf(viewer, "%d component ", Nc);CHKERRQ(ierr);} 157 if (prob->implicit[f]) {ierr = PetscViewerASCIIPrintf(viewer, " (implicit)");CHKERRQ(ierr);} 158 else {ierr = PetscViewerASCIIPrintf(viewer, " (explicit)");CHKERRQ(ierr);} 159 if (prob->adjacency[f*2+0]) { 160 if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM++)");CHKERRQ(ierr);} 161 else {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM)");CHKERRQ(ierr);} 162 } else { 163 if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FEM)");CHKERRQ(ierr);} 164 else {ierr = PetscViewerASCIIPrintf(viewer, " (adj FUNKY)");CHKERRQ(ierr);} 165 } 166 ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 167 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 168 if (id == PETSCFE_CLASSID) {ierr = PetscFEView((PetscFE) obj, viewer);CHKERRQ(ierr);} 169 else if (id == PETSCFV_CLASSID) {ierr = PetscFVView((PetscFV) obj, viewer);CHKERRQ(ierr);} 170 } 171 } 172 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 173 PetscFunctionReturn(0); 174 } 175 176 #undef __FUNCT__ 177 #define __FUNCT__ "PetscDSView" 178 /*@C 179 PetscDSView - Views a PetscDS 180 181 Collective on PetscDS 182 183 Input Parameter: 184 + prob - the PetscDS object to view 185 - v - the viewer 186 187 Level: developer 188 189 .seealso PetscDSDestroy() 190 @*/ 191 PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v) 192 { 193 PetscBool iascii; 194 PetscErrorCode ierr; 195 196 PetscFunctionBegin; 197 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 198 if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);} 199 else {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);} 200 ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 201 if (iascii) {ierr = PetscDSView_Ascii(prob, v);CHKERRQ(ierr);} 202 if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);} 203 PetscFunctionReturn(0); 204 } 205 206 #undef __FUNCT__ 207 #define __FUNCT__ "PetscDSViewFromOptions" 208 /* 209 PetscDSViewFromOptions - Processes command line options to determine if/how a PetscDS is to be viewed. 210 211 Collective on PetscDS 212 213 Input Parameters: 214 + prob - the PetscDS 215 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd' 216 - optionname - option to activate viewing 217 218 Level: intermediate 219 220 .keywords: PetscDS, view, options, database 221 .seealso: VecViewFromOptions(), MatViewFromOptions() 222 */ 223 PetscErrorCode PetscDSViewFromOptions(PetscDS prob, const char prefix[], const char optionname[]) 224 { 225 PetscViewer viewer; 226 PetscViewerFormat format; 227 PetscBool flg; 228 PetscErrorCode ierr; 229 230 PetscFunctionBegin; 231 if (prefix) {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);} 232 else {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), ((PetscObject) prob)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);} 233 if (flg) { 234 ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 235 ierr = PetscDSView(prob, viewer);CHKERRQ(ierr); 236 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 237 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 238 } 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "PetscDSSetFromOptions" 244 /*@ 245 PetscDSSetFromOptions - sets parameters in a PetscDS from the options database 246 247 Collective on PetscDS 248 249 Input Parameter: 250 . prob - the PetscDS object to set options for 251 252 Options Database: 253 254 Level: developer 255 256 .seealso PetscDSView() 257 @*/ 258 PetscErrorCode PetscDSSetFromOptions(PetscDS prob) 259 { 260 const char *defaultType; 261 char name[256]; 262 PetscBool flg; 263 PetscErrorCode ierr; 264 265 PetscFunctionBegin; 266 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 267 if (!((PetscObject) prob)->type_name) { 268 defaultType = PETSCDSBASIC; 269 } else { 270 defaultType = ((PetscObject) prob)->type_name; 271 } 272 ierr = PetscDSRegisterAll();CHKERRQ(ierr); 273 274 ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr); 275 ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr); 276 if (flg) { 277 ierr = PetscDSSetType(prob, name);CHKERRQ(ierr); 278 } else if (!((PetscObject) prob)->type_name) { 279 ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr); 280 } 281 if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);} 282 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 283 ierr = PetscObjectProcessOptionsHandlers((PetscObject) prob);CHKERRQ(ierr); 284 ierr = PetscOptionsEnd();CHKERRQ(ierr); 285 ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr); 286 PetscFunctionReturn(0); 287 } 288 289 #undef __FUNCT__ 290 #define __FUNCT__ "PetscDSSetUp" 291 /*@C 292 PetscDSSetUp - Construct data structures for the PetscDS 293 294 Collective on PetscDS 295 296 Input Parameter: 297 . prob - the PetscDS object to setup 298 299 Level: developer 300 301 .seealso PetscDSView(), PetscDSDestroy() 302 @*/ 303 PetscErrorCode PetscDSSetUp(PetscDS prob) 304 { 305 const PetscInt Nf = prob->Nf; 306 PetscInt dim, work, NcMax = 0, NqMax = 0, f; 307 PetscErrorCode ierr; 308 309 PetscFunctionBegin; 310 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 311 if (prob->setup) PetscFunctionReturn(0); 312 /* Calculate sizes */ 313 ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr); 314 prob->totDim = prob->totDimBd = prob->totComp = 0; 315 ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisBd,Nf,&prob->basisDerBd);CHKERRQ(ierr); 316 for (f = 0; f < Nf; ++f) { 317 PetscFE feBd = (PetscFE) prob->discBd[f]; 318 PetscObject obj; 319 PetscClassId id; 320 PetscQuadrature q; 321 PetscInt Nq = 0, Nb, Nc; 322 323 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 324 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 325 if (id == PETSCFE_CLASSID) { 326 PetscFE fe = (PetscFE) obj; 327 328 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 329 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 330 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 331 ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr); 332 } else if (id == PETSCFV_CLASSID) { 333 PetscFV fv = (PetscFV) obj; 334 335 ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr); 336 Nb = 1; 337 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 338 ierr = PetscFVGetDefaultTabulation(fv, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr); 339 } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 340 if (q) {ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);} 341 NqMax = PetscMax(NqMax, Nq); 342 NcMax = PetscMax(NcMax, Nc); 343 prob->totDim += Nb*Nc; 344 prob->totComp += Nc; 345 if (feBd) { 346 ierr = PetscFEGetDimension(feBd, &Nb);CHKERRQ(ierr); 347 ierr = PetscFEGetNumComponents(feBd, &Nc);CHKERRQ(ierr); 348 ierr = PetscFEGetDefaultTabulation(feBd, &prob->basisBd[f], &prob->basisDerBd[f], NULL);CHKERRQ(ierr); 349 prob->totDimBd += Nb*Nc; 350 } 351 } 352 work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim)); 353 /* Allocate works space */ 354 ierr = PetscMalloc5(prob->totComp,&prob->u,prob->totComp,&prob->u_t,prob->totComp*dim,&prob->u_x,dim,&prob->x,work,&prob->refSpaceDer);CHKERRQ(ierr); 355 ierr = PetscMalloc6(NqMax*NcMax,&prob->f0,NqMax*NcMax*dim,&prob->f1,NqMax*NcMax*NcMax,&prob->g0,NqMax*NcMax*NcMax*dim,&prob->g1,NqMax*NcMax*NcMax*dim,&prob->g2,NqMax*NcMax*NcMax*dim*dim,&prob->g3);CHKERRQ(ierr); 356 if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);} 357 prob->setup = PETSC_TRUE; 358 PetscFunctionReturn(0); 359 } 360 361 #undef __FUNCT__ 362 #define __FUNCT__ "PetscDSDestroyStructs_Static" 363 static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob) 364 { 365 PetscErrorCode ierr; 366 367 PetscFunctionBegin; 368 ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisBd,prob->basisDerBd);CHKERRQ(ierr); 369 ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr); 370 ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr); 371 PetscFunctionReturn(0); 372 } 373 374 #undef __FUNCT__ 375 #define __FUNCT__ "PetscDSEnlarge_Static" 376 static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew) 377 { 378 PetscObject *tmpd, *tmpdbd; 379 PetscBool *tmpi, *tmpa; 380 PointFunc *tmpobj, *tmpf, *tmpg; 381 BdPointFunc *tmpfbd, *tmpgbd; 382 RiemannFunc *tmpr; 383 void **tmpctx; 384 PetscInt Nf = prob->Nf, f, i; 385 PetscErrorCode ierr; 386 387 PetscFunctionBegin; 388 if (Nf >= NfNew) PetscFunctionReturn(0); 389 prob->setup = PETSC_FALSE; 390 ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr); 391 ierr = PetscMalloc4(NfNew, &tmpd, NfNew, &tmpdbd, NfNew, &tmpi, NfNew*2, &tmpa);CHKERRQ(ierr); 392 for (f = 0; f < Nf; ++f) {tmpd[f] = prob->disc[f]; tmpdbd[f] = prob->discBd[f]; tmpi[f] = prob->implicit[f]; for (i = 0; i < 2; ++i) tmpa[f*2+i] = prob->adjacency[f*2+i];} 393 for (f = Nf; f < NfNew; ++f) {ierr = PetscContainerCreate(PetscObjectComm((PetscObject) prob), (PetscContainer *) &tmpd[f]);CHKERRQ(ierr); 394 tmpdbd[f] = NULL; tmpi[f] = PETSC_TRUE; tmpa[f*2+0] = PETSC_FALSE; tmpa[f*2+1] = PETSC_TRUE;} 395 ierr = PetscFree4(prob->disc, prob->discBd, prob->implicit, prob->adjacency);CHKERRQ(ierr); 396 prob->Nf = NfNew; 397 prob->disc = tmpd; 398 prob->discBd = tmpdbd; 399 prob->implicit = tmpi; 400 prob->adjacency = tmpa; 401 ierr = PetscCalloc5(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr); 402 for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f]; 403 for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f]; 404 for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f]; 405 for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f]; 406 for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f]; 407 for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL; 408 for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL; 409 for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL; 410 for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL; 411 for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL; 412 ierr = PetscFree5(prob->obj, prob->f, prob->g, prob->r, prob->ctx);CHKERRQ(ierr); 413 prob->obj = tmpobj; 414 prob->f = tmpf; 415 prob->g = tmpg; 416 prob->r = tmpr; 417 prob->ctx = tmpctx; 418 ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr); 419 for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f]; 420 for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f]; 421 for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL; 422 for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL; 423 ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr); 424 prob->fBd = tmpfbd; 425 prob->gBd = tmpgbd; 426 PetscFunctionReturn(0); 427 } 428 429 #undef __FUNCT__ 430 #define __FUNCT__ "PetscDSDestroy" 431 /*@ 432 PetscDSDestroy - Destroys a PetscDS object 433 434 Collective on PetscDS 435 436 Input Parameter: 437 . prob - the PetscDS object to destroy 438 439 Level: developer 440 441 .seealso PetscDSView() 442 @*/ 443 PetscErrorCode PetscDSDestroy(PetscDS *prob) 444 { 445 PetscInt f; 446 PetscErrorCode ierr; 447 448 PetscFunctionBegin; 449 if (!*prob) PetscFunctionReturn(0); 450 PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1); 451 452 if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);} 453 ((PetscObject) (*prob))->refct = 0; 454 ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr); 455 for (f = 0; f < (*prob)->Nf; ++f) { 456 ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr); 457 ierr = PetscObjectDereference((*prob)->discBd[f]);CHKERRQ(ierr); 458 } 459 ierr = PetscFree4((*prob)->disc, (*prob)->discBd, (*prob)->implicit, (*prob)->adjacency);CHKERRQ(ierr); 460 ierr = PetscFree5((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr); 461 ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr); 462 if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);} 463 ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr); 464 PetscFunctionReturn(0); 465 } 466 467 #undef __FUNCT__ 468 #define __FUNCT__ "PetscDSCreate" 469 /*@ 470 PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType(). 471 472 Collective on MPI_Comm 473 474 Input Parameter: 475 . comm - The communicator for the PetscDS object 476 477 Output Parameter: 478 . prob - The PetscDS object 479 480 Level: beginner 481 482 .seealso: PetscDSSetType(), PETSCDSBASIC 483 @*/ 484 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob) 485 { 486 PetscDS p; 487 PetscErrorCode ierr; 488 489 PetscFunctionBegin; 490 PetscValidPointer(prob, 2); 491 *prob = NULL; 492 ierr = PetscDSInitializePackage();CHKERRQ(ierr); 493 494 ierr = PetscHeaderCreate(p, _p_PetscDS, struct _PetscDSOps, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr); 495 ierr = PetscMemzero(p->ops, sizeof(struct _PetscDSOps));CHKERRQ(ierr); 496 497 p->Nf = 0; 498 p->setup = PETSC_FALSE; 499 500 *prob = p; 501 PetscFunctionReturn(0); 502 } 503 504 #undef __FUNCT__ 505 #define __FUNCT__ "PetscDSGetNumFields" 506 /*@ 507 PetscDSGetNumFields - Returns the number of fields in the DS 508 509 Not collective 510 511 Input Parameter: 512 . prob - The PetscDS object 513 514 Output Parameter: 515 . Nf - The number of fields 516 517 Level: beginner 518 519 .seealso: PetscDSGetSpatialDimension(), PetscDSCreate() 520 @*/ 521 PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf) 522 { 523 PetscFunctionBegin; 524 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 525 PetscValidPointer(Nf, 2); 526 *Nf = prob->Nf; 527 PetscFunctionReturn(0); 528 } 529 530 #undef __FUNCT__ 531 #define __FUNCT__ "PetscDSGetSpatialDimension" 532 /*@ 533 PetscDSGetSpatialDimension - Returns the spatial dimension of the DS 534 535 Not collective 536 537 Input Parameter: 538 . prob - The PetscDS object 539 540 Output Parameter: 541 . dim - The spatial dimension 542 543 Level: beginner 544 545 .seealso: PetscDSGetNumFields(), PetscDSCreate() 546 @*/ 547 PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim) 548 { 549 PetscErrorCode ierr; 550 551 PetscFunctionBegin; 552 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 553 PetscValidPointer(dim, 2); 554 *dim = 0; 555 if (prob->Nf) { 556 PetscObject obj; 557 PetscClassId id; 558 559 ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr); 560 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 561 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);} 562 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);} 563 else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0); 564 } 565 PetscFunctionReturn(0); 566 } 567 568 #undef __FUNCT__ 569 #define __FUNCT__ "PetscDSGetTotalDimension" 570 /*@ 571 PetscDSGetTotalDimension - Returns the total size of the approximation space for this system 572 573 Not collective 574 575 Input Parameter: 576 . prob - The PetscDS object 577 578 Output Parameter: 579 . dim - The total problem dimension 580 581 Level: beginner 582 583 .seealso: PetscDSGetNumFields(), PetscDSCreate() 584 @*/ 585 PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim) 586 { 587 PetscErrorCode ierr; 588 589 PetscFunctionBegin; 590 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 591 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 592 PetscValidPointer(dim, 2); 593 *dim = prob->totDim; 594 PetscFunctionReturn(0); 595 } 596 597 #undef __FUNCT__ 598 #define __FUNCT__ "PetscDSGetTotalBdDimension" 599 /*@ 600 PetscDSGetTotalBdDimension - Returns the total size of the boundary approximation space for this system 601 602 Not collective 603 604 Input Parameter: 605 . prob - The PetscDS object 606 607 Output Parameter: 608 . dim - The total boundary problem dimension 609 610 Level: beginner 611 612 .seealso: PetscDSGetNumFields(), PetscDSCreate() 613 @*/ 614 PetscErrorCode PetscDSGetTotalBdDimension(PetscDS prob, PetscInt *dim) 615 { 616 PetscErrorCode ierr; 617 618 PetscFunctionBegin; 619 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 620 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 621 PetscValidPointer(dim, 2); 622 *dim = prob->totDimBd; 623 PetscFunctionReturn(0); 624 } 625 626 #undef __FUNCT__ 627 #define __FUNCT__ "PetscDSGetTotalComponents" 628 /*@ 629 PetscDSGetTotalComponents - Returns the total number of components in this system 630 631 Not collective 632 633 Input Parameter: 634 . prob - The PetscDS object 635 636 Output Parameter: 637 . dim - The total number of components 638 639 Level: beginner 640 641 .seealso: PetscDSGetNumFields(), PetscDSCreate() 642 @*/ 643 PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc) 644 { 645 PetscErrorCode ierr; 646 647 PetscFunctionBegin; 648 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 649 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 650 PetscValidPointer(Nc, 2); 651 *Nc = prob->totComp; 652 PetscFunctionReturn(0); 653 } 654 655 #undef __FUNCT__ 656 #define __FUNCT__ "PetscDSGetDiscretization" 657 /*@ 658 PetscDSGetDiscretization - Returns the discretization object for the given field 659 660 Not collective 661 662 Input Parameters: 663 + prob - The PetscDS object 664 - f - The field number 665 666 Output Parameter: 667 . disc - The discretization object 668 669 Level: beginner 670 671 .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 672 @*/ 673 PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc) 674 { 675 PetscFunctionBegin; 676 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 677 PetscValidPointer(disc, 3); 678 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 679 *disc = prob->disc[f]; 680 PetscFunctionReturn(0); 681 } 682 683 #undef __FUNCT__ 684 #define __FUNCT__ "PetscDSGetBdDiscretization" 685 /*@ 686 PetscDSGetBdDiscretization - Returns the boundary discretization object for the given field 687 688 Not collective 689 690 Input Parameters: 691 + prob - The PetscDS object 692 - f - The field number 693 694 Output Parameter: 695 . disc - The boundary discretization object 696 697 Level: beginner 698 699 .seealso: PetscDSSetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 700 @*/ 701 PetscErrorCode PetscDSGetBdDiscretization(PetscDS prob, PetscInt f, PetscObject *disc) 702 { 703 PetscFunctionBegin; 704 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 705 PetscValidPointer(disc, 3); 706 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 707 *disc = prob->discBd[f]; 708 PetscFunctionReturn(0); 709 } 710 711 #undef __FUNCT__ 712 #define __FUNCT__ "PetscDSSetDiscretization" 713 /*@ 714 PetscDSSetDiscretization - Sets the discretization object for the given field 715 716 Not collective 717 718 Input Parameters: 719 + prob - The PetscDS object 720 . f - The field number 721 - disc - The discretization object 722 723 Level: beginner 724 725 .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 726 @*/ 727 PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc) 728 { 729 PetscErrorCode ierr; 730 731 PetscFunctionBegin; 732 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 733 PetscValidPointer(disc, 3); 734 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 735 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 736 if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);} 737 prob->disc[f] = disc; 738 ierr = PetscObjectReference(disc);CHKERRQ(ierr); 739 { 740 PetscClassId id; 741 742 ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr); 743 if (id == PETSCFV_CLASSID) { 744 prob->implicit[f] = PETSC_FALSE; 745 prob->adjacency[f*2+0] = PETSC_TRUE; 746 prob->adjacency[f*2+1] = PETSC_FALSE; 747 } 748 } 749 PetscFunctionReturn(0); 750 } 751 752 #undef __FUNCT__ 753 #define __FUNCT__ "PetscDSSetBdDiscretization" 754 /*@ 755 PetscDSSetBdDiscretization - Sets the boundary discretization object for the given field 756 757 Not collective 758 759 Input Parameters: 760 + prob - The PetscDS object 761 . f - The field number 762 - disc - The boundary discretization object 763 764 Level: beginner 765 766 .seealso: PetscDSGetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 767 @*/ 768 PetscErrorCode PetscDSSetBdDiscretization(PetscDS prob, PetscInt f, PetscObject disc) 769 { 770 PetscErrorCode ierr; 771 772 PetscFunctionBegin; 773 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 774 if (disc) PetscValidPointer(disc, 3); 775 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 776 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 777 if (prob->discBd[f]) {ierr = PetscObjectDereference(prob->discBd[f]);CHKERRQ(ierr);} 778 prob->discBd[f] = disc; 779 ierr = PetscObjectReference(disc);CHKERRQ(ierr); 780 PetscFunctionReturn(0); 781 } 782 783 #undef __FUNCT__ 784 #define __FUNCT__ "PetscDSAddDiscretization" 785 /*@ 786 PetscDSAddDiscretization - Adds a discretization object 787 788 Not collective 789 790 Input Parameters: 791 + prob - The PetscDS object 792 - disc - The boundary discretization object 793 794 Level: beginner 795 796 .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 797 @*/ 798 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc) 799 { 800 PetscErrorCode ierr; 801 802 PetscFunctionBegin; 803 ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr); 804 PetscFunctionReturn(0); 805 } 806 807 #undef __FUNCT__ 808 #define __FUNCT__ "PetscDSAddBdDiscretization" 809 /*@ 810 PetscDSAddBdDiscretization - Adds a boundary discretization object 811 812 Not collective 813 814 Input Parameters: 815 + prob - The PetscDS object 816 - disc - The boundary discretization object 817 818 Level: beginner 819 820 .seealso: PetscDSGetBdDiscretization(), PetscDSSetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 821 @*/ 822 PetscErrorCode PetscDSAddBdDiscretization(PetscDS prob, PetscObject disc) 823 { 824 PetscErrorCode ierr; 825 826 PetscFunctionBegin; 827 ierr = PetscDSSetBdDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr); 828 PetscFunctionReturn(0); 829 } 830 831 #undef __FUNCT__ 832 #define __FUNCT__ "PetscDSGetImplicit" 833 /*@ 834 PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX 835 836 Not collective 837 838 Input Parameters: 839 + prob - The PetscDS object 840 - f - The field number 841 842 Output Parameter: 843 . implicit - The flag indicating what kind of solve to use for this field 844 845 Level: developer 846 847 .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 848 @*/ 849 PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit) 850 { 851 PetscFunctionBegin; 852 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 853 PetscValidPointer(implicit, 3); 854 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 855 *implicit = prob->implicit[f]; 856 PetscFunctionReturn(0); 857 } 858 859 #undef __FUNCT__ 860 #define __FUNCT__ "PetscDSSetImplicit" 861 /*@ 862 PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX 863 864 Not collective 865 866 Input Parameters: 867 + prob - The PetscDS object 868 . f - The field number 869 - implicit - The flag indicating what kind of solve to use for this field 870 871 Level: developer 872 873 .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 874 @*/ 875 PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit) 876 { 877 PetscFunctionBegin; 878 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 879 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 880 prob->implicit[f] = implicit; 881 PetscFunctionReturn(0); 882 } 883 884 #undef __FUNCT__ 885 #define __FUNCT__ "PetscDSGetAdjacency" 886 /*@ 887 PetscDSGetAdjacency - Returns the flags for determining variable influence 888 889 Not collective 890 891 Input Parameters: 892 + prob - The PetscDS object 893 - f - The field number 894 895 Output Parameter: 896 + useCone - Flag for variable influence starting with the cone operation 897 - useClosure - Flag for variable influence using transitive closure 898 899 Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure() 900 901 Level: developer 902 903 .seealso: PetscDSSetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 904 @*/ 905 PetscErrorCode PetscDSGetAdjacency(PetscDS prob, PetscInt f, PetscBool *useCone, PetscBool *useClosure) 906 { 907 PetscFunctionBegin; 908 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 909 PetscValidPointer(useCone, 3); 910 PetscValidPointer(useClosure, 4); 911 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 912 *useCone = prob->adjacency[f*2+0]; 913 *useClosure = prob->adjacency[f*2+1]; 914 PetscFunctionReturn(0); 915 } 916 917 #undef __FUNCT__ 918 #define __FUNCT__ "PetscDSSetAdjacency" 919 /*@ 920 PetscDSSetAdjacency - Set the flags for determining variable influence 921 922 Not collective 923 924 Input Parameters: 925 + prob - The PetscDS object 926 . f - The field number 927 . useCone - Flag for variable influence starting with the cone operation 928 - useClosure - Flag for variable influence using transitive closure 929 930 Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure() 931 932 Level: developer 933 934 .seealso: PetscDSGetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 935 @*/ 936 PetscErrorCode PetscDSSetAdjacency(PetscDS prob, PetscInt f, PetscBool useCone, PetscBool useClosure) 937 { 938 PetscFunctionBegin; 939 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 940 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 941 prob->adjacency[f*2+0] = useCone; 942 prob->adjacency[f*2+1] = useClosure; 943 PetscFunctionReturn(0); 944 } 945 946 #undef __FUNCT__ 947 #define __FUNCT__ "PetscDSGetObjective" 948 PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f, 949 void (**obj)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar obj[])) 950 { 951 PetscFunctionBegin; 952 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 953 PetscValidPointer(obj, 2); 954 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 955 *obj = prob->obj[f]; 956 PetscFunctionReturn(0); 957 } 958 959 #undef __FUNCT__ 960 #define __FUNCT__ "PetscDSSetObjective" 961 PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f, 962 void (*obj)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar obj[])) 963 { 964 PetscErrorCode ierr; 965 966 PetscFunctionBegin; 967 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 968 PetscValidFunction(obj, 2); 969 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 970 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 971 prob->obj[f] = obj; 972 PetscFunctionReturn(0); 973 } 974 975 #undef __FUNCT__ 976 #define __FUNCT__ "PetscDSGetResidual" 977 PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f, 978 void (**f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f0[]), 979 void (**f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f1[])) 980 { 981 PetscFunctionBegin; 982 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 983 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 984 if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];} 985 if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];} 986 PetscFunctionReturn(0); 987 } 988 989 #undef __FUNCT__ 990 #define __FUNCT__ "PetscDSSetResidual" 991 PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f, 992 void (*f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f0[]), 993 void (*f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f1[])) 994 { 995 PetscErrorCode ierr; 996 997 PetscFunctionBegin; 998 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 999 PetscValidFunction(f0, 3); 1000 PetscValidFunction(f1, 4); 1001 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1002 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 1003 prob->f[f*2+0] = f0; 1004 prob->f[f*2+1] = f1; 1005 PetscFunctionReturn(0); 1006 } 1007 1008 #undef __FUNCT__ 1009 #define __FUNCT__ "PetscDSGetJacobian" 1010 PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g, 1011 void (**g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g0[]), 1012 void (**g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g1[]), 1013 void (**g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g2[]), 1014 void (**g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g3[])) 1015 { 1016 PetscFunctionBegin; 1017 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1018 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1019 if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf); 1020 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];} 1021 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];} 1022 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];} 1023 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];} 1024 PetscFunctionReturn(0); 1025 } 1026 1027 #undef __FUNCT__ 1028 #define __FUNCT__ "PetscDSSetJacobian" 1029 PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g, 1030 void (*g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g0[]), 1031 void (*g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g1[]), 1032 void (*g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g2[]), 1033 void (*g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g3[])) 1034 { 1035 PetscErrorCode ierr; 1036 1037 PetscFunctionBegin; 1038 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1039 if (g0) PetscValidFunction(g0, 4); 1040 if (g1) PetscValidFunction(g1, 5); 1041 if (g2) PetscValidFunction(g2, 6); 1042 if (g3) PetscValidFunction(g3, 7); 1043 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1044 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 1045 ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr); 1046 prob->g[(f*prob->Nf + g)*4+0] = g0; 1047 prob->g[(f*prob->Nf + g)*4+1] = g1; 1048 prob->g[(f*prob->Nf + g)*4+2] = g2; 1049 prob->g[(f*prob->Nf + g)*4+3] = g3; 1050 PetscFunctionReturn(0); 1051 } 1052 1053 #undef __FUNCT__ 1054 #define __FUNCT__ "PetscDSGetRiemannSolver" 1055 /*@C 1056 PetscDSGetRiemannSolver - Returns the Riemann solver for the given field 1057 1058 Not collective 1059 1060 Input Arguments: 1061 + prob - The PetscDS object 1062 - f - The field number 1063 1064 Output Argument: 1065 . r - Riemann solver 1066 1067 Calling sequence for r: 1068 1069 $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 1070 1071 + x - The coordinates at a point on the interface 1072 . n - The normal vector to the interface 1073 . uL - The state vector to the left of the interface 1074 . uR - The state vector to the right of the interface 1075 . flux - output array of flux through the interface 1076 - ctx - optional user context 1077 1078 Level: intermediate 1079 1080 .seealso: PetscDSSetRiemannSolver() 1081 @*/ 1082 PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f, 1083 void (**r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)) 1084 { 1085 PetscFunctionBegin; 1086 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1087 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1088 PetscValidPointer(r, 3); 1089 *r = prob->r[f]; 1090 PetscFunctionReturn(0); 1091 } 1092 1093 #undef __FUNCT__ 1094 #define __FUNCT__ "PetscDSSetRiemannSolver" 1095 /*@C 1096 PetscDSSetRiemannSolver - Sets the Riemann solver for the given field 1097 1098 Not collective 1099 1100 Input Arguments: 1101 + prob - The PetscDS object 1102 . f - The field number 1103 - r - Riemann solver 1104 1105 Calling sequence for r: 1106 1107 $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 1108 1109 + x - The coordinates at a point on the interface 1110 . n - The normal vector to the interface 1111 . uL - The state vector to the left of the interface 1112 . uR - The state vector to the right of the interface 1113 . flux - output array of flux through the interface 1114 - ctx - optional user context 1115 1116 Level: intermediate 1117 1118 .seealso: PetscDSGetRiemannSolver() 1119 @*/ 1120 PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f, 1121 void (*r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)) 1122 { 1123 PetscErrorCode ierr; 1124 1125 PetscFunctionBegin; 1126 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1127 PetscValidFunction(r, 3); 1128 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1129 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 1130 prob->r[f] = r; 1131 PetscFunctionReturn(0); 1132 } 1133 1134 #undef __FUNCT__ 1135 #define __FUNCT__ "PetscDSGetContext" 1136 PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx) 1137 { 1138 PetscFunctionBegin; 1139 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1140 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1141 PetscValidPointer(ctx, 3); 1142 *ctx = prob->ctx[f]; 1143 PetscFunctionReturn(0); 1144 } 1145 1146 #undef __FUNCT__ 1147 #define __FUNCT__ "PetscDSSetContext" 1148 PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx) 1149 { 1150 PetscErrorCode ierr; 1151 1152 PetscFunctionBegin; 1153 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1154 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1155 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 1156 prob->ctx[f] = ctx; 1157 PetscFunctionReturn(0); 1158 } 1159 1160 #undef __FUNCT__ 1161 #define __FUNCT__ "PetscDSGetBdResidual" 1162 PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f, 1163 void (**f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f0[]), 1164 void (**f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f1[])) 1165 { 1166 PetscFunctionBegin; 1167 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1168 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1169 if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];} 1170 if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];} 1171 PetscFunctionReturn(0); 1172 } 1173 1174 #undef __FUNCT__ 1175 #define __FUNCT__ "PetscDSSetBdResidual" 1176 PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f, 1177 void (*f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f0[]), 1178 void (*f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f1[])) 1179 { 1180 PetscErrorCode ierr; 1181 1182 PetscFunctionBegin; 1183 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1184 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1185 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 1186 if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;} 1187 if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;} 1188 PetscFunctionReturn(0); 1189 } 1190 1191 #undef __FUNCT__ 1192 #define __FUNCT__ "PetscDSGetBdJacobian" 1193 PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g, 1194 void (**g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g0[]), 1195 void (**g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g1[]), 1196 void (**g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g2[]), 1197 void (**g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g3[])) 1198 { 1199 PetscFunctionBegin; 1200 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1201 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1202 if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf); 1203 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];} 1204 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];} 1205 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];} 1206 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];} 1207 PetscFunctionReturn(0); 1208 } 1209 1210 #undef __FUNCT__ 1211 #define __FUNCT__ "PetscDSSetBdJacobian" 1212 PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g, 1213 void (*g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g0[]), 1214 void (*g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g1[]), 1215 void (*g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g2[]), 1216 void (*g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g3[])) 1217 { 1218 PetscErrorCode ierr; 1219 1220 PetscFunctionBegin; 1221 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1222 if (g0) PetscValidFunction(g0, 4); 1223 if (g1) PetscValidFunction(g1, 5); 1224 if (g2) PetscValidFunction(g2, 6); 1225 if (g3) PetscValidFunction(g3, 7); 1226 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1227 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 1228 ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr); 1229 prob->gBd[(f*prob->Nf + g)*4+0] = g0; 1230 prob->gBd[(f*prob->Nf + g)*4+1] = g1; 1231 prob->gBd[(f*prob->Nf + g)*4+2] = g2; 1232 prob->gBd[(f*prob->Nf + g)*4+3] = g3; 1233 PetscFunctionReturn(0); 1234 } 1235 1236 #undef __FUNCT__ 1237 #define __FUNCT__ "PetscDSGetFieldOffset" 1238 /*@ 1239 PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis 1240 1241 Not collective 1242 1243 Input Parameters: 1244 + prob - The PetscDS object 1245 - f - The field number 1246 1247 Output Parameter: 1248 . off - The offset 1249 1250 Level: beginner 1251 1252 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate() 1253 @*/ 1254 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off) 1255 { 1256 PetscInt g; 1257 PetscErrorCode ierr; 1258 1259 PetscFunctionBegin; 1260 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1261 PetscValidPointer(off, 3); 1262 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1263 *off = 0; 1264 for (g = 0; g < f; ++g) { 1265 PetscFE fe = (PetscFE) prob->disc[g]; 1266 PetscInt Nb, Nc; 1267 1268 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1269 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1270 *off += Nb*Nc; 1271 } 1272 PetscFunctionReturn(0); 1273 } 1274 1275 #undef __FUNCT__ 1276 #define __FUNCT__ "PetscDSGetBdFieldOffset" 1277 /*@ 1278 PetscDSGetBdFieldOffset - Returns the offset of the given field in the full space boundary basis 1279 1280 Not collective 1281 1282 Input Parameters: 1283 + prob - The PetscDS object 1284 - f - The field number 1285 1286 Output Parameter: 1287 . off - The boundary offset 1288 1289 Level: beginner 1290 1291 .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate() 1292 @*/ 1293 PetscErrorCode PetscDSGetBdFieldOffset(PetscDS prob, PetscInt f, PetscInt *off) 1294 { 1295 PetscInt g; 1296 PetscErrorCode ierr; 1297 1298 PetscFunctionBegin; 1299 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1300 PetscValidPointer(off, 3); 1301 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1302 *off = 0; 1303 for (g = 0; g < f; ++g) { 1304 PetscFE fe = (PetscFE) prob->discBd[g]; 1305 PetscInt Nb, Nc; 1306 1307 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1308 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1309 *off += Nb*Nc; 1310 } 1311 PetscFunctionReturn(0); 1312 } 1313 1314 #undef __FUNCT__ 1315 #define __FUNCT__ "PetscDSGetComponentOffset" 1316 /*@ 1317 PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point 1318 1319 Not collective 1320 1321 Input Parameters: 1322 + prob - The PetscDS object 1323 - f - The field number 1324 1325 Output Parameter: 1326 . off - The offset 1327 1328 Level: beginner 1329 1330 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate() 1331 @*/ 1332 PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off) 1333 { 1334 PetscInt g; 1335 PetscErrorCode ierr; 1336 1337 PetscFunctionBegin; 1338 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1339 PetscValidPointer(off, 3); 1340 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1341 *off = 0; 1342 for (g = 0; g < f; ++g) { 1343 PetscFE fe = (PetscFE) prob->disc[g]; 1344 PetscInt Nc; 1345 1346 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1347 *off += Nc; 1348 } 1349 PetscFunctionReturn(0); 1350 } 1351 1352 #undef __FUNCT__ 1353 #define __FUNCT__ "PetscDSGetTabulation" 1354 /*@C 1355 PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization 1356 1357 Not collective 1358 1359 Input Parameter: 1360 . prob - The PetscDS object 1361 1362 Output Parameters: 1363 + basis - The basis function tabulation at quadrature points 1364 - basisDer - The basis function derivative tabulation at quadrature points 1365 1366 Level: intermediate 1367 1368 .seealso: PetscDSGetBdTabulation(), PetscDSCreate() 1369 @*/ 1370 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer) 1371 { 1372 PetscErrorCode ierr; 1373 1374 PetscFunctionBegin; 1375 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1376 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 1377 if (basis) {PetscValidPointer(basis, 2); *basis = prob->basis;} 1378 if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;} 1379 PetscFunctionReturn(0); 1380 } 1381 1382 #undef __FUNCT__ 1383 #define __FUNCT__ "PetscDSGetBdTabulation" 1384 /*@C 1385 PetscDSGetBdTabulation - Return the basis tabulation at quadrature points for the boundary discretization 1386 1387 Not collective 1388 1389 Input Parameter: 1390 . prob - The PetscDS object 1391 1392 Output Parameters: 1393 + basis - The basis function tabulation at quadrature points 1394 - basisDer - The basis function derivative tabulation at quadrature points 1395 1396 Level: intermediate 1397 1398 .seealso: PetscDSGetTabulation(), PetscDSCreate() 1399 @*/ 1400 PetscErrorCode PetscDSGetBdTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer) 1401 { 1402 PetscErrorCode ierr; 1403 1404 PetscFunctionBegin; 1405 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1406 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 1407 if (basis) {PetscValidPointer(basis, 2); *basis = prob->basisBd;} 1408 if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerBd;} 1409 PetscFunctionReturn(0); 1410 } 1411 1412 #undef __FUNCT__ 1413 #define __FUNCT__ "PetscDSGetEvaluationArrays" 1414 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x) 1415 { 1416 PetscErrorCode ierr; 1417 1418 PetscFunctionBegin; 1419 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1420 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 1421 if (u) {PetscValidPointer(u, 2); *u = prob->u;} 1422 if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;} 1423 if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;} 1424 PetscFunctionReturn(0); 1425 } 1426 1427 #undef __FUNCT__ 1428 #define __FUNCT__ "PetscDSGetWeakFormArrays" 1429 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3) 1430 { 1431 PetscErrorCode ierr; 1432 1433 PetscFunctionBegin; 1434 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1435 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 1436 if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;} 1437 if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;} 1438 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;} 1439 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;} 1440 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;} 1441 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;} 1442 PetscFunctionReturn(0); 1443 } 1444 1445 #undef __FUNCT__ 1446 #define __FUNCT__ "PetscDSGetRefCoordArrays" 1447 PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer) 1448 { 1449 PetscErrorCode ierr; 1450 1451 PetscFunctionBegin; 1452 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1453 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 1454 if (x) {PetscValidPointer(x, 2); *x = prob->x;} 1455 if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;} 1456 PetscFunctionReturn(0); 1457 } 1458 1459 #undef __FUNCT__ 1460 #define __FUNCT__ "PetscDSDestroy_Basic" 1461 static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob) 1462 { 1463 PetscFunctionBegin; 1464 PetscFunctionReturn(0); 1465 } 1466 1467 #undef __FUNCT__ 1468 #define __FUNCT__ "PetscDSInitialize_Basic" 1469 static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob) 1470 { 1471 PetscFunctionBegin; 1472 prob->ops->setfromoptions = NULL; 1473 prob->ops->setup = NULL; 1474 prob->ops->view = NULL; 1475 prob->ops->destroy = PetscDSDestroy_Basic; 1476 PetscFunctionReturn(0); 1477 } 1478 1479 /*MC 1480 PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions 1481 1482 Level: intermediate 1483 1484 .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType() 1485 M*/ 1486 1487 #undef __FUNCT__ 1488 #define __FUNCT__ "PetscDSCreate_Basic" 1489 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob) 1490 { 1491 PetscDS_Basic *b; 1492 PetscErrorCode ierr; 1493 1494 PetscFunctionBegin; 1495 PetscValidHeaderSpecific(prob, PETSCSPACE_CLASSID, 1); 1496 ierr = PetscNewLog(prob, &b);CHKERRQ(ierr); 1497 prob->data = b; 1498 1499 ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr); 1500 PetscFunctionReturn(0); 1501 } 1502