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 /* A PetscDS (Discrete System) encodes a set of equations posed in a discrete space, which represents a set of 9 nonlinear continuum equations. The equations can have multiple fields, each field having a different 10 discretization. In addition, different pieces of the domain can have different field combinations and equations. 11 12 The DS provides the user a description of the approximation space on any given cell. It also gives pointwise 13 functions representing the equations. 14 15 Each field is associated with a label, marking the cells on which it is supported. Note that a field can be 16 supported on the closure of a cell not in the label due to overlap of the boundary of neighboring cells. The DM 17 then creates a DS for each set of cells with identical approximation spaces. When assembling, the user asks for 18 the space associated with a given cell. DMPlex uses the labels associated with each DS in the default integration loop. 19 */ 20 21 /*@C 22 PetscDSRegister - Adds a new `PetscDS` implementation 23 24 Not Collective 25 26 Input Parameters: 27 + name - The name of a new user-defined creation routine 28 - create_func - The creation routine itself 29 30 Sample usage: 31 .vb 32 PetscDSRegister("my_ds", MyPetscDSCreate); 33 .ve 34 35 Then, your PetscDS type can be chosen with the procedural interface via 36 .vb 37 PetscDSCreate(MPI_Comm, PetscDS *); 38 PetscDSSetType(PetscDS, "my_ds"); 39 .ve 40 or at runtime via the option 41 .vb 42 -petscds_type my_ds 43 .ve 44 45 Level: advanced 46 47 Note: 48 `PetscDSRegister()` may be called multiple times to add several user-defined `PetscDSs` 49 50 Fortran Note: 51 Not available from Fortran 52 53 .seealso: `PetscDSType`, `PetscDS`, `PetscDSRegisterAll()`, `PetscDSRegisterDestroy()` 54 @*/ 55 PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS)) 56 { 57 PetscFunctionBegin; 58 PetscCall(PetscFunctionListAdd(&PetscDSList, sname, function)); 59 PetscFunctionReturn(PETSC_SUCCESS); 60 } 61 62 /*@C 63 PetscDSSetType - Builds a particular `PetscDS` 64 65 Collective on prob 66 67 Input Parameters: 68 + prob - The `PetscDS` object 69 - name - The `PetscDSType` 70 71 Options Database Key: 72 . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types 73 74 Level: intermediate 75 76 Fortran Note: 77 Not available from Fortran 78 79 .seealso: `PetscDSType`, `PetscDS`, `PetscDSGetType()`, `PetscDSCreate()` 80 @*/ 81 PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name) 82 { 83 PetscErrorCode (*r)(PetscDS); 84 PetscBool match; 85 86 PetscFunctionBegin; 87 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 88 PetscCall(PetscObjectTypeCompare((PetscObject)prob, name, &match)); 89 if (match) PetscFunctionReturn(PETSC_SUCCESS); 90 91 PetscCall(PetscDSRegisterAll()); 92 PetscCall(PetscFunctionListFind(PetscDSList, name, &r)); 93 PetscCheck(r, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name); 94 95 PetscTryTypeMethod(prob, destroy); 96 prob->ops->destroy = NULL; 97 98 PetscCall((*r)(prob)); 99 PetscCall(PetscObjectChangeTypeName((PetscObject)prob, name)); 100 PetscFunctionReturn(PETSC_SUCCESS); 101 } 102 103 /*@C 104 PetscDSGetType - Gets the `PetscDSType` name (as a string) from the `PetscDS` 105 106 Not Collective 107 108 Input Parameter: 109 . prob - The `PetscDS` 110 111 Output Parameter: 112 . name - The `PetscDSType` name 113 114 Level: intermediate 115 116 Fortran Note: 117 Not available from Fortran 118 119 .seealso: `PetscDSType`, `PetscDS`, `PetscDSSetType()`, `PetscDSCreate()` 120 @*/ 121 PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name) 122 { 123 PetscFunctionBegin; 124 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 125 PetscValidPointer(name, 2); 126 PetscCall(PetscDSRegisterAll()); 127 *name = ((PetscObject)prob)->type_name; 128 PetscFunctionReturn(PETSC_SUCCESS); 129 } 130 131 static PetscErrorCode PetscDSView_Ascii(PetscDS ds, PetscViewer viewer) 132 { 133 PetscViewerFormat format; 134 const PetscScalar *constants; 135 PetscInt Nf, numConstants, f; 136 137 PetscFunctionBegin; 138 PetscCall(PetscDSGetNumFields(ds, &Nf)); 139 PetscCall(PetscViewerGetFormat(viewer, &format)); 140 PetscCall(PetscViewerASCIIPrintf(viewer, "Discrete System with %" PetscInt_FMT " fields\n", Nf)); 141 PetscCall(PetscViewerASCIIPushTab(viewer)); 142 PetscCall(PetscViewerASCIIPrintf(viewer, " cell total dim %" PetscInt_FMT " total comp %" PetscInt_FMT "\n", ds->totDim, ds->totComp)); 143 if (ds->isCohesive) PetscCall(PetscViewerASCIIPrintf(viewer, " cohesive cell\n")); 144 for (f = 0; f < Nf; ++f) { 145 DSBoundary b; 146 PetscObject obj; 147 PetscClassId id; 148 PetscQuadrature q; 149 const char *name; 150 PetscInt Nc, Nq, Nqc; 151 152 PetscCall(PetscDSGetDiscretization(ds, f, &obj)); 153 PetscCall(PetscObjectGetClassId(obj, &id)); 154 PetscCall(PetscObjectGetName(obj, &name)); 155 PetscCall(PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>")); 156 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 157 if (id == PETSCFE_CLASSID) { 158 PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 159 PetscCall(PetscFEGetQuadrature((PetscFE)obj, &q)); 160 PetscCall(PetscViewerASCIIPrintf(viewer, " FEM")); 161 } else if (id == PETSCFV_CLASSID) { 162 PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc)); 163 PetscCall(PetscFVGetQuadrature((PetscFV)obj, &q)); 164 PetscCall(PetscViewerASCIIPrintf(viewer, " FVM")); 165 } else SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f); 166 if (Nc > 1) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " components", Nc)); 167 else PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " component ", Nc)); 168 if (ds->implicit[f]) PetscCall(PetscViewerASCIIPrintf(viewer, " (implicit)")); 169 else PetscCall(PetscViewerASCIIPrintf(viewer, " (explicit)")); 170 if (q) { 171 PetscCall(PetscQuadratureGetData(q, NULL, &Nqc, &Nq, NULL, NULL)); 172 PetscCall(PetscViewerASCIIPrintf(viewer, " (Nq %" PetscInt_FMT " Nqc %" PetscInt_FMT ")", Nq, Nqc)); 173 } 174 PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT "-jet", ds->jetDegree[f])); 175 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 176 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 177 PetscCall(PetscViewerASCIIPushTab(viewer)); 178 if (id == PETSCFE_CLASSID) PetscCall(PetscFEView((PetscFE)obj, viewer)); 179 else if (id == PETSCFV_CLASSID) PetscCall(PetscFVView((PetscFV)obj, viewer)); 180 PetscCall(PetscViewerASCIIPopTab(viewer)); 181 182 for (b = ds->boundary; b; b = b->next) { 183 char *name; 184 PetscInt c, i; 185 186 if (b->field != f) continue; 187 PetscCall(PetscViewerASCIIPushTab(viewer)); 188 PetscCall(PetscViewerASCIIPrintf(viewer, "Boundary %s (%s) %s\n", b->name, b->lname, DMBoundaryConditionTypes[b->type])); 189 if (!b->Nc) { 190 PetscCall(PetscViewerASCIIPrintf(viewer, " all components\n")); 191 } else { 192 PetscCall(PetscViewerASCIIPrintf(viewer, " components: ")); 193 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 194 for (c = 0; c < b->Nc; ++c) { 195 if (c > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 196 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT, b->comps[c])); 197 } 198 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 199 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 200 } 201 PetscCall(PetscViewerASCIIPrintf(viewer, " values: ")); 202 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 203 for (i = 0; i < b->Nv; ++i) { 204 if (i > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 205 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT, b->values[i])); 206 } 207 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 208 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 209 #if defined(__clang__) 210 PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wformat-pedantic"); 211 #elif defined(__GNUC__) || defined(__GNUG__) 212 PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wformat"); 213 #endif 214 if (b->func) { 215 PetscCall(PetscDLAddr(b->func, &name)); 216 if (name) PetscCall(PetscViewerASCIIPrintf(viewer, " func: %s\n", name)); 217 else PetscCall(PetscViewerASCIIPrintf(viewer, " func: %p\n", b->func)); 218 PetscCall(PetscFree(name)); 219 } 220 if (b->func_t) { 221 PetscCall(PetscDLAddr(b->func_t, &name)); 222 if (name) PetscCall(PetscViewerASCIIPrintf(viewer, " func_t: %s\n", name)); 223 else PetscCall(PetscViewerASCIIPrintf(viewer, " func_t: %p\n", b->func_t)); 224 PetscCall(PetscFree(name)); 225 } 226 PETSC_PRAGMA_DIAGNOSTIC_IGNORED_END(); 227 PetscCall(PetscWeakFormView(b->wf, viewer)); 228 PetscCall(PetscViewerASCIIPopTab(viewer)); 229 } 230 } 231 PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 232 if (numConstants) { 233 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " constants\n", numConstants)); 234 PetscCall(PetscViewerASCIIPushTab(viewer)); 235 for (f = 0; f < numConstants; ++f) PetscCall(PetscViewerASCIIPrintf(viewer, "%g\n", (double)PetscRealPart(constants[f]))); 236 PetscCall(PetscViewerASCIIPopTab(viewer)); 237 } 238 PetscCall(PetscWeakFormView(ds->wf, viewer)); 239 PetscCall(PetscViewerASCIIPopTab(viewer)); 240 PetscFunctionReturn(PETSC_SUCCESS); 241 } 242 243 /*@C 244 PetscDSViewFromOptions - View a `PetscDS` based on values in the options database 245 246 Collective on PetscDS 247 248 Input Parameters: 249 + A - the `PetscDS` object 250 . obj - Optional object 251 - name - command line option 252 253 Level: intermediate 254 255 .seealso: `PetscDSType`, `PetscDS`, `PetscDSView()`, `PetscObjectViewFromOptions()`, `PetscDSCreate()` 256 @*/ 257 PetscErrorCode PetscDSViewFromOptions(PetscDS A, PetscObject obj, const char name[]) 258 { 259 PetscFunctionBegin; 260 PetscValidHeaderSpecific(A, PETSCDS_CLASSID, 1); 261 PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 262 PetscFunctionReturn(PETSC_SUCCESS); 263 } 264 265 /*@C 266 PetscDSView - Views a `PetscDS` 267 268 Collective on prob 269 270 Input Parameters: 271 + prob - the `PetscDS` object to view 272 - v - the viewer 273 274 Level: developer 275 276 .seealso: `PetscDSType`, `PetscDS`, `PetscViewer`, `PetscDSDestroy()` 277 @*/ 278 PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v) 279 { 280 PetscBool iascii; 281 282 PetscFunctionBegin; 283 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 284 if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)prob), &v)); 285 else PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2); 286 PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii)); 287 if (iascii) PetscCall(PetscDSView_Ascii(prob, v)); 288 PetscTryTypeMethod(prob, view, v); 289 PetscFunctionReturn(PETSC_SUCCESS); 290 } 291 292 /*@ 293 PetscDSSetFromOptions - sets parameters in a `PetscDS` from the options database 294 295 Collective on prob 296 297 Input Parameter: 298 . prob - the `PetscDS` object to set options for 299 300 Options Database Keys: 301 + -petscds_type <type> - Set the `PetscDS` type 302 . -petscds_view <view opt> - View the `PetscDS` 303 . -petscds_jac_pre - Turn formation of a separate Jacobian preconditioner on or off 304 . -bc_<name> <ids> - Specify a list of label ids for a boundary condition 305 - -bc_<name>_comp <comps> - Specify a list of field components to constrain for a boundary condition 306 307 Level: intermediate 308 309 .seealso: `PetscDS`, `PetscDSView()` 310 @*/ 311 PetscErrorCode PetscDSSetFromOptions(PetscDS prob) 312 { 313 DSBoundary b; 314 const char *defaultType; 315 char name[256]; 316 PetscBool flg; 317 318 PetscFunctionBegin; 319 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 320 if (!((PetscObject)prob)->type_name) { 321 defaultType = PETSCDSBASIC; 322 } else { 323 defaultType = ((PetscObject)prob)->type_name; 324 } 325 PetscCall(PetscDSRegisterAll()); 326 327 PetscObjectOptionsBegin((PetscObject)prob); 328 for (b = prob->boundary; b; b = b->next) { 329 char optname[1024]; 330 PetscInt ids[1024], len = 1024; 331 PetscBool flg; 332 333 PetscCall(PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name)); 334 PetscCall(PetscMemzero(ids, sizeof(ids))); 335 PetscCall(PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg)); 336 if (flg) { 337 b->Nv = len; 338 PetscCall(PetscFree(b->values)); 339 PetscCall(PetscMalloc1(len, &b->values)); 340 PetscCall(PetscArraycpy(b->values, ids, len)); 341 PetscCall(PetscWeakFormRewriteKeys(b->wf, b->label, len, b->values)); 342 } 343 len = 1024; 344 PetscCall(PetscSNPrintf(optname, sizeof(optname), "-bc_%s_comp", b->name)); 345 PetscCall(PetscMemzero(ids, sizeof(ids))); 346 PetscCall(PetscOptionsIntArray(optname, "List of boundary field components", "", ids, &len, &flg)); 347 if (flg) { 348 b->Nc = len; 349 PetscCall(PetscFree(b->comps)); 350 PetscCall(PetscMalloc1(len, &b->comps)); 351 PetscCall(PetscArraycpy(b->comps, ids, len)); 352 } 353 } 354 PetscCall(PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg)); 355 if (flg) { 356 PetscCall(PetscDSSetType(prob, name)); 357 } else if (!((PetscObject)prob)->type_name) { 358 PetscCall(PetscDSSetType(prob, defaultType)); 359 } 360 PetscCall(PetscOptionsBool("-petscds_jac_pre", "Discrete System", "PetscDSUseJacobianPreconditioner", prob->useJacPre, &prob->useJacPre, &flg)); 361 PetscCall(PetscOptionsBool("-petscds_force_quad", "Discrete System", "PetscDSSetForceQuad", prob->forceQuad, &prob->forceQuad, &flg)); 362 PetscTryTypeMethod(prob, setfromoptions); 363 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 364 PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)prob, PetscOptionsObject)); 365 PetscOptionsEnd(); 366 if (prob->Nf) PetscCall(PetscDSViewFromOptions(prob, NULL, "-petscds_view")); 367 PetscFunctionReturn(PETSC_SUCCESS); 368 } 369 370 /*@C 371 PetscDSSetUp - Construct data structures for the `PetscDS` 372 373 Collective on prob 374 375 Input Parameter: 376 . prob - the `PetscDS` object to setup 377 378 Level: developer 379 380 .seealso: `PetscDS`, `PetscDSView()`, `PetscDSDestroy()` 381 @*/ 382 PetscErrorCode PetscDSSetUp(PetscDS prob) 383 { 384 const PetscInt Nf = prob->Nf; 385 PetscBool hasH = PETSC_FALSE; 386 PetscInt dim, dimEmbed, NbMax = 0, NcMax = 0, NqMax = 0, NsMax = 1, f; 387 PetscInt gorder = -1, fgorder = -1; 388 389 PetscFunctionBegin; 390 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 391 if (prob->setup) PetscFunctionReturn(PETSC_SUCCESS); 392 /* Calculate sizes */ 393 PetscCall(PetscDSGetSpatialDimension(prob, &dim)); 394 PetscCall(PetscDSGetCoordinateDimension(prob, &dimEmbed)); 395 prob->totDim = prob->totComp = 0; 396 PetscCall(PetscMalloc2(Nf, &prob->Nc, Nf, &prob->Nb)); 397 PetscCall(PetscCalloc2(Nf + 1, &prob->off, Nf + 1, &prob->offDer)); 398 PetscCall(PetscCalloc6(Nf + 1, &prob->offCohesive[0], Nf + 1, &prob->offCohesive[1], Nf + 1, &prob->offCohesive[2], Nf + 1, &prob->offDerCohesive[0], Nf + 1, &prob->offDerCohesive[1], Nf + 1, &prob->offDerCohesive[2])); 399 PetscCall(PetscMalloc2(Nf, &prob->T, Nf, &prob->Tf)); 400 if (prob->forceQuad) { 401 PetscQuadrature mq = NULL, mfq = NULL; 402 PetscInt maxOrder = -1, maxFOrder = -1; 403 404 for (f = 0; f < Nf; ++f) { 405 PetscObject obj; 406 PetscClassId id; 407 PetscQuadrature q = NULL, fq = NULL; 408 PetscInt order = -1, forder = -1; 409 410 PetscCall(PetscDSGetDiscretization(prob, f, &obj)); 411 if (!obj) continue; 412 PetscCall(PetscObjectGetClassId(obj, &id)); 413 if (id == PETSCFE_CLASSID) { 414 PetscFE fe = (PetscFE)obj; 415 416 PetscCall(PetscFEGetQuadrature(fe, &q)); 417 PetscCall(PetscFEGetFaceQuadrature(fe, &fq)); 418 } else if (id == PETSCFV_CLASSID) { 419 PetscFV fv = (PetscFV)obj; 420 421 PetscCall(PetscFVGetQuadrature(fv, &q)); 422 } 423 if (q) PetscCall(PetscQuadratureGetOrder(q, &order)); 424 if (fq) PetscCall(PetscQuadratureGetOrder(fq, &forder)); 425 if (order > maxOrder) { 426 maxOrder = order; 427 mq = q; 428 } 429 if (forder > maxFOrder) { 430 maxFOrder = forder; 431 mfq = fq; 432 } 433 } 434 for (f = 0; f < Nf; ++f) { 435 PetscObject obj; 436 PetscClassId id; 437 438 PetscCall(PetscDSGetDiscretization(prob, f, &obj)); 439 if (!obj) continue; 440 PetscCall(PetscObjectGetClassId(obj, &id)); 441 if (id == PETSCFE_CLASSID) { 442 PetscFE fe = (PetscFE)obj; 443 444 PetscCall(PetscFESetQuadrature(fe, mq)); 445 PetscCall(PetscFESetFaceQuadrature(fe, mfq)); 446 } else if (id == PETSCFV_CLASSID) { 447 PetscFV fv = (PetscFV)obj; 448 449 PetscCall(PetscFVSetQuadrature(fv, mq)); 450 } 451 } 452 } 453 for (f = 0; f < Nf; ++f) { 454 PetscObject obj; 455 PetscClassId id; 456 PetscQuadrature q = NULL; 457 PetscInt Nq = 0, Nb, Nc; 458 459 PetscCall(PetscDSGetDiscretization(prob, f, &obj)); 460 if (prob->jetDegree[f] > 1) hasH = PETSC_TRUE; 461 if (!obj) { 462 /* Empty mesh */ 463 Nb = Nc = 0; 464 prob->T[f] = prob->Tf[f] = NULL; 465 } else { 466 PetscCall(PetscObjectGetClassId(obj, &id)); 467 if (id == PETSCFE_CLASSID) { 468 PetscFE fe = (PetscFE)obj; 469 470 PetscCall(PetscFEGetQuadrature(fe, &q)); 471 { 472 PetscQuadrature fq; 473 PetscInt order, forder; 474 475 PetscCall(PetscQuadratureGetOrder(q, &order)); 476 if (gorder < 0) gorder = order; 477 PetscCheck(order == gorder, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field %" PetscInt_FMT " cell quadrature order %" PetscInt_FMT " != %" PetscInt_FMT " DS cell quadrature order", f, order, gorder); 478 PetscCall(PetscFEGetFaceQuadrature(fe, &fq)); 479 if (fq) { 480 PetscCall(PetscQuadratureGetOrder(fq, &forder)); 481 if (fgorder < 0) fgorder = forder; 482 PetscCheck(forder == fgorder, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field %" PetscInt_FMT " face quadrature order %" PetscInt_FMT " != %" PetscInt_FMT " DS face quadrature order", f, forder, fgorder); 483 } 484 } 485 PetscCall(PetscFEGetDimension(fe, &Nb)); 486 PetscCall(PetscFEGetNumComponents(fe, &Nc)); 487 PetscCall(PetscFEGetCellTabulation(fe, prob->jetDegree[f], &prob->T[f])); 488 PetscCall(PetscFEGetFaceTabulation(fe, prob->jetDegree[f], &prob->Tf[f])); 489 } else if (id == PETSCFV_CLASSID) { 490 PetscFV fv = (PetscFV)obj; 491 492 PetscCall(PetscFVGetQuadrature(fv, &q)); 493 PetscCall(PetscFVGetNumComponents(fv, &Nc)); 494 Nb = Nc; 495 PetscCall(PetscFVGetCellTabulation(fv, &prob->T[f])); 496 /* TODO: should PetscFV also have face tabulation? Otherwise there will be a null pointer in prob->basisFace */ 497 } else SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f); 498 } 499 prob->Nc[f] = Nc; 500 prob->Nb[f] = Nb; 501 prob->off[f + 1] = Nc + prob->off[f]; 502 prob->offDer[f + 1] = Nc * dim + prob->offDer[f]; 503 prob->offCohesive[0][f + 1] = (prob->cohesive[f] ? Nc : Nc * 2) + prob->offCohesive[0][f]; 504 prob->offDerCohesive[0][f + 1] = (prob->cohesive[f] ? Nc : Nc * 2) * dimEmbed + prob->offDerCohesive[0][f]; 505 prob->offCohesive[1][f] = (prob->cohesive[f] ? 0 : Nc) + prob->offCohesive[0][f]; 506 prob->offDerCohesive[1][f] = (prob->cohesive[f] ? 0 : Nc) * dimEmbed + prob->offDerCohesive[0][f]; 507 prob->offCohesive[2][f + 1] = (prob->cohesive[f] ? Nc : Nc * 2) + prob->offCohesive[2][f]; 508 prob->offDerCohesive[2][f + 1] = (prob->cohesive[f] ? Nc : Nc * 2) * dimEmbed + prob->offDerCohesive[2][f]; 509 if (q) PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL)); 510 NqMax = PetscMax(NqMax, Nq); 511 NbMax = PetscMax(NbMax, Nb); 512 NcMax = PetscMax(NcMax, Nc); 513 prob->totDim += Nb; 514 prob->totComp += Nc; 515 /* There are two faces for all fields on a cohesive cell, except for cohesive fields */ 516 if (prob->isCohesive && !prob->cohesive[f]) prob->totDim += Nb; 517 } 518 prob->offCohesive[1][Nf] = prob->offCohesive[0][Nf]; 519 prob->offDerCohesive[1][Nf] = prob->offDerCohesive[0][Nf]; 520 /* Allocate works space */ 521 NsMax = 2; /* A non-cohesive discretizations can be used on a cohesive cell, so we need this extra workspace for all DS */ 522 PetscCall(PetscMalloc3(NsMax * prob->totComp, &prob->u, NsMax * prob->totComp, &prob->u_t, NsMax * prob->totComp * dimEmbed + (hasH ? NsMax * prob->totComp * dimEmbed * dimEmbed : 0), &prob->u_x)); 523 PetscCall(PetscMalloc5(dimEmbed, &prob->x, NbMax * NcMax, &prob->basisReal, NbMax * NcMax * dimEmbed, &prob->basisDerReal, NbMax * NcMax, &prob->testReal, NbMax * NcMax * dimEmbed, &prob->testDerReal)); 524 PetscCall(PetscMalloc6(NsMax * NqMax * NcMax, &prob->f0, NsMax * NqMax * NcMax * dimEmbed, &prob->f1, NsMax * NsMax * NqMax * NcMax * NcMax, &prob->g0, NsMax * NsMax * NqMax * NcMax * NcMax * dimEmbed, &prob->g1, NsMax * NsMax * NqMax * NcMax * NcMax * dimEmbed, 525 &prob->g2, NsMax * NsMax * NqMax * NcMax * NcMax * dimEmbed * dimEmbed, &prob->g3)); 526 PetscTryTypeMethod(prob, setup); 527 prob->setup = PETSC_TRUE; 528 PetscFunctionReturn(PETSC_SUCCESS); 529 } 530 531 static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob) 532 { 533 PetscFunctionBegin; 534 PetscCall(PetscFree2(prob->Nc, prob->Nb)); 535 PetscCall(PetscFree2(prob->off, prob->offDer)); 536 PetscCall(PetscFree6(prob->offCohesive[0], prob->offCohesive[1], prob->offCohesive[2], prob->offDerCohesive[0], prob->offDerCohesive[1], prob->offDerCohesive[2])); 537 PetscCall(PetscFree2(prob->T, prob->Tf)); 538 PetscCall(PetscFree3(prob->u, prob->u_t, prob->u_x)); 539 PetscCall(PetscFree5(prob->x, prob->basisReal, prob->basisDerReal, prob->testReal, prob->testDerReal)); 540 PetscCall(PetscFree6(prob->f0, prob->f1, prob->g0, prob->g1, prob->g2, prob->g3)); 541 PetscFunctionReturn(PETSC_SUCCESS); 542 } 543 544 static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew) 545 { 546 PetscObject *tmpd; 547 PetscBool *tmpi; 548 PetscInt *tmpk; 549 PetscBool *tmpc; 550 PetscPointFunc *tmpup; 551 PetscSimplePointFunc *tmpexactSol, *tmpexactSol_t; 552 void **tmpexactCtx, **tmpexactCtx_t; 553 void **tmpctx; 554 PetscInt Nf = prob->Nf, f; 555 556 PetscFunctionBegin; 557 if (Nf >= NfNew) PetscFunctionReturn(PETSC_SUCCESS); 558 prob->setup = PETSC_FALSE; 559 PetscCall(PetscDSDestroyStructs_Static(prob)); 560 PetscCall(PetscMalloc4(NfNew, &tmpd, NfNew, &tmpi, NfNew, &tmpc, NfNew, &tmpk)); 561 for (f = 0; f < Nf; ++f) { 562 tmpd[f] = prob->disc[f]; 563 tmpi[f] = prob->implicit[f]; 564 tmpc[f] = prob->cohesive[f]; 565 tmpk[f] = prob->jetDegree[f]; 566 } 567 for (f = Nf; f < NfNew; ++f) { 568 tmpd[f] = NULL; 569 tmpi[f] = PETSC_TRUE, tmpc[f] = PETSC_FALSE; 570 tmpk[f] = 1; 571 } 572 PetscCall(PetscFree4(prob->disc, prob->implicit, prob->cohesive, prob->jetDegree)); 573 PetscCall(PetscWeakFormSetNumFields(prob->wf, NfNew)); 574 prob->Nf = NfNew; 575 prob->disc = tmpd; 576 prob->implicit = tmpi; 577 prob->cohesive = tmpc; 578 prob->jetDegree = tmpk; 579 PetscCall(PetscCalloc2(NfNew, &tmpup, NfNew, &tmpctx)); 580 for (f = 0; f < Nf; ++f) tmpup[f] = prob->update[f]; 581 for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f]; 582 for (f = Nf; f < NfNew; ++f) tmpup[f] = NULL; 583 for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL; 584 PetscCall(PetscFree2(prob->update, prob->ctx)); 585 prob->update = tmpup; 586 prob->ctx = tmpctx; 587 PetscCall(PetscCalloc4(NfNew, &tmpexactSol, NfNew, &tmpexactCtx, NfNew, &tmpexactSol_t, NfNew, &tmpexactCtx_t)); 588 for (f = 0; f < Nf; ++f) tmpexactSol[f] = prob->exactSol[f]; 589 for (f = 0; f < Nf; ++f) tmpexactCtx[f] = prob->exactCtx[f]; 590 for (f = 0; f < Nf; ++f) tmpexactSol_t[f] = prob->exactSol_t[f]; 591 for (f = 0; f < Nf; ++f) tmpexactCtx_t[f] = prob->exactCtx_t[f]; 592 for (f = Nf; f < NfNew; ++f) tmpexactSol[f] = NULL; 593 for (f = Nf; f < NfNew; ++f) tmpexactCtx[f] = NULL; 594 for (f = Nf; f < NfNew; ++f) tmpexactSol_t[f] = NULL; 595 for (f = Nf; f < NfNew; ++f) tmpexactCtx_t[f] = NULL; 596 PetscCall(PetscFree4(prob->exactSol, prob->exactCtx, prob->exactSol_t, prob->exactCtx_t)); 597 prob->exactSol = tmpexactSol; 598 prob->exactCtx = tmpexactCtx; 599 prob->exactSol_t = tmpexactSol_t; 600 prob->exactCtx_t = tmpexactCtx_t; 601 PetscFunctionReturn(PETSC_SUCCESS); 602 } 603 604 /*@ 605 PetscDSDestroy - Destroys a PetscDS object 606 607 Collective on prob 608 609 Input Parameter: 610 . prob - the PetscDS object to destroy 611 612 Level: developer 613 614 .seealso: `PetscDSView()` 615 @*/ 616 PetscErrorCode PetscDSDestroy(PetscDS *ds) 617 { 618 PetscInt f; 619 620 PetscFunctionBegin; 621 if (!*ds) PetscFunctionReturn(PETSC_SUCCESS); 622 PetscValidHeaderSpecific((*ds), PETSCDS_CLASSID, 1); 623 624 if (--((PetscObject)(*ds))->refct > 0) { 625 *ds = NULL; 626 PetscFunctionReturn(PETSC_SUCCESS); 627 } 628 ((PetscObject)(*ds))->refct = 0; 629 if ((*ds)->subprobs) { 630 PetscInt dim, d; 631 632 PetscCall(PetscDSGetSpatialDimension(*ds, &dim)); 633 for (d = 0; d < dim; ++d) PetscCall(PetscDSDestroy(&(*ds)->subprobs[d])); 634 } 635 PetscCall(PetscFree((*ds)->subprobs)); 636 PetscCall(PetscDSDestroyStructs_Static(*ds)); 637 for (f = 0; f < (*ds)->Nf; ++f) PetscCall(PetscObjectDereference((*ds)->disc[f])); 638 PetscCall(PetscFree4((*ds)->disc, (*ds)->implicit, (*ds)->cohesive, (*ds)->jetDegree)); 639 PetscCall(PetscWeakFormDestroy(&(*ds)->wf)); 640 PetscCall(PetscFree2((*ds)->update, (*ds)->ctx)); 641 PetscCall(PetscFree4((*ds)->exactSol, (*ds)->exactCtx, (*ds)->exactSol_t, (*ds)->exactCtx_t)); 642 PetscTryTypeMethod((*ds), destroy); 643 PetscCall(PetscDSDestroyBoundary(*ds)); 644 PetscCall(PetscFree((*ds)->constants)); 645 PetscCall(PetscHeaderDestroy(ds)); 646 PetscFunctionReturn(PETSC_SUCCESS); 647 } 648 649 /*@ 650 PetscDSCreate - Creates an empty `PetscDS` object. The type can then be set with `PetscDSSetType()`. 651 652 Collective 653 654 Input Parameter: 655 . comm - The communicator for the `PetscDS` object 656 657 Output Parameter: 658 . ds - The `PetscDS` object 659 660 Level: beginner 661 662 .seealso: `PetscDS`, `PetscDSSetType()`, `PETSCDSBASIC`, `PetscDSType` 663 @*/ 664 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *ds) 665 { 666 PetscDS p; 667 668 PetscFunctionBegin; 669 PetscValidPointer(ds, 2); 670 *ds = NULL; 671 PetscCall(PetscDSInitializePackage()); 672 673 PetscCall(PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView)); 674 675 p->Nf = 0; 676 p->setup = PETSC_FALSE; 677 p->numConstants = 0; 678 p->constants = NULL; 679 p->dimEmbed = -1; 680 p->useJacPre = PETSC_TRUE; 681 p->forceQuad = PETSC_TRUE; 682 PetscCall(PetscWeakFormCreate(comm, &p->wf)); 683 684 *ds = p; 685 PetscFunctionReturn(PETSC_SUCCESS); 686 } 687 688 /*@ 689 PetscDSGetNumFields - Returns the number of fields in the `PetscDS` 690 691 Not collective 692 693 Input Parameter: 694 . prob - The PetscDS object 695 696 Output Parameter: 697 . Nf - The number of fields 698 699 Level: beginner 700 701 .seealso: `PetscDS`, `PetscDSGetSpatialDimension()`, `PetscDSCreate()` 702 @*/ 703 PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf) 704 { 705 PetscFunctionBegin; 706 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 707 PetscValidIntPointer(Nf, 2); 708 *Nf = prob->Nf; 709 PetscFunctionReturn(PETSC_SUCCESS); 710 } 711 712 /*@ 713 PetscDSGetSpatialDimension - Returns the spatial dimension of the `PetscDS`, meaning the topological dimension of the discretizations 714 715 Not collective 716 717 Input Parameter: 718 . prob - The `PetscDS` object 719 720 Output Parameter: 721 . dim - The spatial dimension 722 723 Level: beginner 724 725 .seealso: `PetscDS`, `PetscDSGetCoordinateDimension()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 726 @*/ 727 PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim) 728 { 729 PetscFunctionBegin; 730 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 731 PetscValidIntPointer(dim, 2); 732 *dim = 0; 733 if (prob->Nf) { 734 PetscObject obj; 735 PetscClassId id; 736 737 PetscCall(PetscDSGetDiscretization(prob, 0, &obj)); 738 if (obj) { 739 PetscCall(PetscObjectGetClassId(obj, &id)); 740 if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetSpatialDimension((PetscFE)obj, dim)); 741 else if (id == PETSCFV_CLASSID) PetscCall(PetscFVGetSpatialDimension((PetscFV)obj, dim)); 742 else SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0); 743 } 744 } 745 PetscFunctionReturn(PETSC_SUCCESS); 746 } 747 748 /*@ 749 PetscDSGetCoordinateDimension - Returns the coordinate dimension of the `PetscDS`, meaning the dimension of the space into which the discretiaztions are embedded 750 751 Not collective 752 753 Input Parameter: 754 . prob - The `PetscDS` object 755 756 Output Parameter: 757 . dimEmbed - The coordinate dimension 758 759 Level: beginner 760 761 .seealso: `PetscDS`, `PetscDSSetCoordinateDimension()`, `PetscDSGetSpatialDimension()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 762 @*/ 763 PetscErrorCode PetscDSGetCoordinateDimension(PetscDS prob, PetscInt *dimEmbed) 764 { 765 PetscFunctionBegin; 766 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 767 PetscValidIntPointer(dimEmbed, 2); 768 PetscCheck(prob->dimEmbed >= 0, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONGSTATE, "No coordinate dimension set for this DS"); 769 *dimEmbed = prob->dimEmbed; 770 PetscFunctionReturn(PETSC_SUCCESS); 771 } 772 773 /*@ 774 PetscDSSetCoordinateDimension - Set the coordinate dimension of the `PetscDS`, meaning the dimension of the space into which the discretiaztions are embedded 775 776 Logically collective on prob 777 778 Input Parameters: 779 + prob - The `PetscDS` object 780 - dimEmbed - The coordinate dimension 781 782 Level: beginner 783 784 .seealso: `PetscDS`, `PetscDSGetCoordinateDimension()`, `PetscDSGetSpatialDimension()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 785 @*/ 786 PetscErrorCode PetscDSSetCoordinateDimension(PetscDS prob, PetscInt dimEmbed) 787 { 788 PetscFunctionBegin; 789 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 790 PetscCheck(dimEmbed >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coordinate dimension must be non-negative, not %" PetscInt_FMT, dimEmbed); 791 prob->dimEmbed = dimEmbed; 792 PetscFunctionReturn(PETSC_SUCCESS); 793 } 794 795 /*@ 796 PetscDSGetForceQuad - Returns the flag to force matching quadratures among the field discretizations 797 798 Not collective 799 800 Input Parameter: 801 . prob - The `PetscDS` object 802 803 Output Parameter: 804 . forceQuad - The flag 805 806 Level: intermediate 807 808 .seealso: `PetscDS`, `PetscDSSetForceQuad()`, `PetscDSGetDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 809 @*/ 810 PetscErrorCode PetscDSGetForceQuad(PetscDS ds, PetscBool *forceQuad) 811 { 812 PetscFunctionBegin; 813 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 814 PetscValidIntPointer(forceQuad, 2); 815 *forceQuad = ds->forceQuad; 816 PetscFunctionReturn(PETSC_SUCCESS); 817 } 818 819 /*@ 820 PetscDSSetForceQuad - Set the flag to force matching quadratures among the field discretizations 821 822 Logically collective on ds 823 824 Input Parameters: 825 + ds - The `PetscDS` object 826 - forceQuad - The flag 827 828 Level: intermediate 829 830 .seealso: `PetscDS`, `PetscDSGetForceQuad()`, `PetscDSGetDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 831 @*/ 832 PetscErrorCode PetscDSSetForceQuad(PetscDS ds, PetscBool forceQuad) 833 { 834 PetscFunctionBegin; 835 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 836 ds->forceQuad = forceQuad; 837 PetscFunctionReturn(PETSC_SUCCESS); 838 } 839 840 /*@ 841 PetscDSIsCohesive - Returns the flag indicating that this `PetscDS` is for a cohesive cell 842 843 Not collective 844 845 Input Parameter: 846 . ds - The `PetscDS` object 847 848 Output Parameter: 849 . isCohesive - The flag 850 851 Level: developer 852 853 .seealso: `PetscDS`, `PetscDSGetNumCohesive()`, `PetscDSGetCohesive()`, `PetscDSSetCohesive()`, `PetscDSCreate()` 854 @*/ 855 PetscErrorCode PetscDSIsCohesive(PetscDS ds, PetscBool *isCohesive) 856 { 857 PetscFunctionBegin; 858 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 859 PetscValidBoolPointer(isCohesive, 2); 860 *isCohesive = ds->isCohesive; 861 PetscFunctionReturn(PETSC_SUCCESS); 862 } 863 864 /*@ 865 PetscDSGetNumCohesive - Returns the numer of cohesive fields, meaning those defined on the interior of a cohesive cell 866 867 Not collective 868 869 Input Parameter: 870 . ds - The `PetscDS` object 871 872 Output Parameter: 873 . numCohesive - The number of cohesive fields 874 875 Level: developer 876 877 .seealso: `PetscDS`, `PetscDSSetCohesive()`, `PetscDSCreate()` 878 @*/ 879 PetscErrorCode PetscDSGetNumCohesive(PetscDS ds, PetscInt *numCohesive) 880 { 881 PetscInt f; 882 883 PetscFunctionBegin; 884 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 885 PetscValidIntPointer(numCohesive, 2); 886 *numCohesive = 0; 887 for (f = 0; f < ds->Nf; ++f) *numCohesive += ds->cohesive[f] ? 1 : 0; 888 PetscFunctionReturn(PETSC_SUCCESS); 889 } 890 891 /*@ 892 PetscDSGetCohesive - Returns the flag indicating that a field is cohesive, meaning it is defined on the interior of a cohesive cell 893 894 Not collective 895 896 Input Parameters: 897 + ds - The `PetscDS` object 898 - f - The field index 899 900 Output Parameter: 901 . isCohesive - The flag 902 903 Level: developer 904 905 .seealso: `PetscDS`, `PetscDSSetCohesive()`, `PetscDSIsCohesive()`, `PetscDSCreate()` 906 @*/ 907 PetscErrorCode PetscDSGetCohesive(PetscDS ds, PetscInt f, PetscBool *isCohesive) 908 { 909 PetscFunctionBegin; 910 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 911 PetscValidBoolPointer(isCohesive, 3); 912 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 913 *isCohesive = ds->cohesive[f]; 914 PetscFunctionReturn(PETSC_SUCCESS); 915 } 916 917 /*@ 918 PetscDSSetCohesive - Set the flag indicating that a field is cohesive, meaning it is defined on the interior of a cohesive cell 919 920 Not collective 921 922 Input Parameters: 923 + ds - The `PetscDS` object 924 . f - The field index 925 - isCohesive - The flag for a cohesive field 926 927 Level: developer 928 929 .seealso: `PetscDS`, `PetscDSGetCohesive()`, `PetscDSIsCohesive()`, `PetscDSCreate()` 930 @*/ 931 PetscErrorCode PetscDSSetCohesive(PetscDS ds, PetscInt f, PetscBool isCohesive) 932 { 933 PetscInt i; 934 935 PetscFunctionBegin; 936 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 937 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 938 ds->cohesive[f] = isCohesive; 939 ds->isCohesive = PETSC_FALSE; 940 for (i = 0; i < ds->Nf; ++i) ds->isCohesive = ds->isCohesive || ds->cohesive[f] ? PETSC_TRUE : PETSC_FALSE; 941 PetscFunctionReturn(PETSC_SUCCESS); 942 } 943 944 /*@ 945 PetscDSGetTotalDimension - Returns the total size of the approximation space for this system 946 947 Not collective 948 949 Input Parameter: 950 . prob - The `PetscDS` object 951 952 Output Parameter: 953 . dim - The total problem dimension 954 955 Level: beginner 956 957 .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()` 958 @*/ 959 PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim) 960 { 961 PetscFunctionBegin; 962 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 963 PetscCall(PetscDSSetUp(prob)); 964 PetscValidIntPointer(dim, 2); 965 *dim = prob->totDim; 966 PetscFunctionReturn(PETSC_SUCCESS); 967 } 968 969 /*@ 970 PetscDSGetTotalComponents - Returns the total number of components in this system 971 972 Not collective 973 974 Input Parameter: 975 . prob - The `PetscDS` object 976 977 Output Parameter: 978 . dim - The total number of components 979 980 Level: beginner 981 982 .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()` 983 @*/ 984 PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc) 985 { 986 PetscFunctionBegin; 987 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 988 PetscCall(PetscDSSetUp(prob)); 989 PetscValidIntPointer(Nc, 2); 990 *Nc = prob->totComp; 991 PetscFunctionReturn(PETSC_SUCCESS); 992 } 993 994 /*@ 995 PetscDSGetDiscretization - Returns the discretization object for the given field 996 997 Not collective 998 999 Input Parameters: 1000 + prob - The `PetscDS` object 1001 - f - The field number 1002 1003 Output Parameter: 1004 . disc - The discretization object 1005 1006 Level: beginner 1007 1008 .seealso: `PetscDS`, `PetscFE`, `PetscFV`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 1009 @*/ 1010 PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc) 1011 { 1012 PetscFunctionBeginHot; 1013 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1014 PetscValidPointer(disc, 3); 1015 PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf); 1016 *disc = prob->disc[f]; 1017 PetscFunctionReturn(PETSC_SUCCESS); 1018 } 1019 1020 /*@ 1021 PetscDSSetDiscretization - Sets the discretization object for the given field 1022 1023 Not collective 1024 1025 Input Parameters: 1026 + prob - The `PetscDS` object 1027 . f - The field number 1028 - disc - The discretization object 1029 1030 Level: beginner 1031 1032 .seealso: `PetscDS`, `PetscFE`, `PetscFV`, `PetscDSGetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 1033 @*/ 1034 PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc) 1035 { 1036 PetscFunctionBegin; 1037 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1038 if (disc) PetscValidPointer(disc, 3); 1039 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 1040 PetscCall(PetscDSEnlarge_Static(prob, f + 1)); 1041 PetscCall(PetscObjectDereference(prob->disc[f])); 1042 prob->disc[f] = disc; 1043 PetscCall(PetscObjectReference(disc)); 1044 if (disc) { 1045 PetscClassId id; 1046 1047 PetscCall(PetscObjectGetClassId(disc, &id)); 1048 if (id == PETSCFE_CLASSID) { 1049 PetscCall(PetscDSSetImplicit(prob, f, PETSC_TRUE)); 1050 } else if (id == PETSCFV_CLASSID) { 1051 PetscCall(PetscDSSetImplicit(prob, f, PETSC_FALSE)); 1052 } 1053 PetscCall(PetscDSSetJetDegree(prob, f, 1)); 1054 } 1055 PetscFunctionReturn(PETSC_SUCCESS); 1056 } 1057 1058 /*@ 1059 PetscDSGetWeakForm - Returns the weak form object 1060 1061 Not collective 1062 1063 Input Parameter: 1064 . ds - The `PetscDS` object 1065 1066 Output Parameter: 1067 . wf - The weak form object 1068 1069 Level: beginner 1070 1071 .seealso: `PetscWeakForm`, `PetscDSSetWeakForm()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 1072 @*/ 1073 PetscErrorCode PetscDSGetWeakForm(PetscDS ds, PetscWeakForm *wf) 1074 { 1075 PetscFunctionBegin; 1076 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1077 PetscValidPointer(wf, 2); 1078 *wf = ds->wf; 1079 PetscFunctionReturn(PETSC_SUCCESS); 1080 } 1081 1082 /*@ 1083 PetscDSSetWeakForm - Sets the weak form object 1084 1085 Not collective 1086 1087 Input Parameters: 1088 + ds - The `PetscDS` object 1089 - wf - The weak form object 1090 1091 Level: beginner 1092 1093 .seealso: `PetscWeakForm`, `PetscDSGetWeakForm()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 1094 @*/ 1095 PetscErrorCode PetscDSSetWeakForm(PetscDS ds, PetscWeakForm wf) 1096 { 1097 PetscFunctionBegin; 1098 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1099 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 2); 1100 PetscCall(PetscObjectDereference((PetscObject)ds->wf)); 1101 ds->wf = wf; 1102 PetscCall(PetscObjectReference((PetscObject)wf)); 1103 PetscCall(PetscWeakFormSetNumFields(wf, ds->Nf)); 1104 PetscFunctionReturn(PETSC_SUCCESS); 1105 } 1106 1107 /*@ 1108 PetscDSAddDiscretization - Adds a discretization object 1109 1110 Not collective 1111 1112 Input Parameters: 1113 + prob - The `PetscDS` object 1114 - disc - The boundary discretization object 1115 1116 Level: beginner 1117 1118 .seealso: `PetscWeakForm`, `PetscDSGetDiscretization()`, `PetscDSSetDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 1119 @*/ 1120 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc) 1121 { 1122 PetscFunctionBegin; 1123 PetscCall(PetscDSSetDiscretization(prob, prob->Nf, disc)); 1124 PetscFunctionReturn(PETSC_SUCCESS); 1125 } 1126 1127 /*@ 1128 PetscDSGetQuadrature - Returns the quadrature, which must agree for all fields in the `PetscDS` 1129 1130 Not collective 1131 1132 Input Parameter: 1133 . prob - The `PetscDS` object 1134 1135 Output Parameter: 1136 . q - The quadrature object 1137 1138 Level: intermediate 1139 1140 .seealso: `PetscDS`, `PetscQuadrature`, `PetscDSSetImplicit()`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 1141 @*/ 1142 PetscErrorCode PetscDSGetQuadrature(PetscDS prob, PetscQuadrature *q) 1143 { 1144 PetscObject obj; 1145 PetscClassId id; 1146 1147 PetscFunctionBegin; 1148 *q = NULL; 1149 if (!prob->Nf) PetscFunctionReturn(PETSC_SUCCESS); 1150 PetscCall(PetscDSGetDiscretization(prob, 0, &obj)); 1151 PetscCall(PetscObjectGetClassId(obj, &id)); 1152 if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetQuadrature((PetscFE)obj, q)); 1153 else if (id == PETSCFV_CLASSID) PetscCall(PetscFVGetQuadrature((PetscFV)obj, q)); 1154 else SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0); 1155 PetscFunctionReturn(PETSC_SUCCESS); 1156 } 1157 1158 /*@ 1159 PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for `TSIMEX` 1160 1161 Not collective 1162 1163 Input Parameters: 1164 + prob - The `PetscDS` object 1165 - f - The field number 1166 1167 Output Parameter: 1168 . implicit - The flag indicating what kind of solve to use for this field 1169 1170 Level: developer 1171 1172 .seealso: `TSIMEX`, `PetscDS`, `PetscDSSetImplicit()`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 1173 @*/ 1174 PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit) 1175 { 1176 PetscFunctionBegin; 1177 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1178 PetscValidBoolPointer(implicit, 3); 1179 PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf); 1180 *implicit = prob->implicit[f]; 1181 PetscFunctionReturn(PETSC_SUCCESS); 1182 } 1183 1184 /*@ 1185 PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for `TSIMEX` 1186 1187 Not collective 1188 1189 Input Parameters: 1190 + prob - The `PetscDS` object 1191 . f - The field number 1192 - implicit - The flag indicating what kind of solve to use for this field 1193 1194 Level: developer 1195 1196 .seealso: `TSIMEX`, `PetscDSGetImplicit()`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 1197 @*/ 1198 PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit) 1199 { 1200 PetscFunctionBegin; 1201 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1202 PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf); 1203 prob->implicit[f] = implicit; 1204 PetscFunctionReturn(PETSC_SUCCESS); 1205 } 1206 1207 /*@ 1208 PetscDSGetJetDegree - Returns the highest derivative for this field equation, or the k-jet that the discretization needs to tabulate. 1209 1210 Not collective 1211 1212 Input Parameters: 1213 + ds - The `PetscDS` object 1214 - f - The field number 1215 1216 Output Parameter: 1217 . k - The highest derivative we need to tabulate 1218 1219 Level: developer 1220 1221 .seealso: `PetscDS`, `PetscDSSetJetDegree()`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 1222 @*/ 1223 PetscErrorCode PetscDSGetJetDegree(PetscDS ds, PetscInt f, PetscInt *k) 1224 { 1225 PetscFunctionBegin; 1226 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1227 PetscValidIntPointer(k, 3); 1228 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 1229 *k = ds->jetDegree[f]; 1230 PetscFunctionReturn(PETSC_SUCCESS); 1231 } 1232 1233 /*@ 1234 PetscDSSetJetDegree - Set the highest derivative for this field equation, or the k-jet that the discretization needs to tabulate. 1235 1236 Not collective 1237 1238 Input Parameters: 1239 + ds - The `PetscDS` object 1240 . f - The field number 1241 - k - The highest derivative we need to tabulate 1242 1243 Level: developer 1244 1245 .seealso: ``PetscDS`, PetscDSGetJetDegree()`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 1246 @*/ 1247 PetscErrorCode PetscDSSetJetDegree(PetscDS ds, PetscInt f, PetscInt k) 1248 { 1249 PetscFunctionBegin; 1250 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1251 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 1252 ds->jetDegree[f] = k; 1253 PetscFunctionReturn(PETSC_SUCCESS); 1254 } 1255 1256 PetscErrorCode PetscDSGetObjective(PetscDS ds, PetscInt f, void (**obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])) 1257 { 1258 PetscPointFunc *tmp; 1259 PetscInt n; 1260 1261 PetscFunctionBegin; 1262 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1263 PetscValidPointer(obj, 3); 1264 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 1265 PetscCall(PetscWeakFormGetObjective(ds->wf, NULL, 0, f, 0, &n, &tmp)); 1266 *obj = tmp ? tmp[0] : NULL; 1267 PetscFunctionReturn(PETSC_SUCCESS); 1268 } 1269 1270 PetscErrorCode PetscDSSetObjective(PetscDS ds, PetscInt f, void (*obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])) 1271 { 1272 PetscFunctionBegin; 1273 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1274 if (obj) PetscValidFunction(obj, 3); 1275 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 1276 PetscCall(PetscWeakFormSetIndexObjective(ds->wf, NULL, 0, f, 0, 0, obj)); 1277 PetscFunctionReturn(PETSC_SUCCESS); 1278 } 1279 1280 /*@C 1281 PetscDSGetResidual - Get the pointwise residual function for a given test field 1282 1283 Not collective 1284 1285 Input Parameters: 1286 + ds - The `PetscDS` 1287 - f - The test field number 1288 1289 Output Parameters: 1290 + f0 - integrand for the test function term 1291 - f1 - integrand for the test function gradient term 1292 1293 Calling sequence for the callbacks f0 and f1: 1294 .vb 1295 f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1296 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1297 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1298 PetscReal t, const PetscReal x[], PetscScalar f0[]) 1299 .ve 1300 + dim - the spatial dimension 1301 . Nf - the number of fields 1302 . uOff - the offset into u[] and u_t[] for each field 1303 . uOff_x - the offset into u_x[] for each field 1304 . u - each field evaluated at the current point 1305 . u_t - the time derivative of each field evaluated at the current point 1306 . u_x - the gradient of each field evaluated at the current point 1307 . aOff - the offset into a[] and a_t[] for each auxiliary field 1308 . aOff_x - the offset into a_x[] for each auxiliary field 1309 . a - each auxiliary field evaluated at the current point 1310 . a_t - the time derivative of each auxiliary field evaluated at the current point 1311 . a_x - the gradient of auxiliary each field evaluated at the current point 1312 . t - current time 1313 . x - coordinates of the current point 1314 . numConstants - number of constant parameters 1315 . constants - constant parameters 1316 - f0 - output values at the current point 1317 1318 Level: intermediate 1319 1320 Note: 1321 We are using a first order FEM model for the weak form: \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t) 1322 1323 .seealso: `PetscDS`, `PetscDSSetResidual()` 1324 @*/ 1325 PetscErrorCode PetscDSGetResidual(PetscDS ds, PetscInt f, void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 1326 { 1327 PetscPointFunc *tmp0, *tmp1; 1328 PetscInt n0, n1; 1329 1330 PetscFunctionBegin; 1331 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1332 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 1333 PetscCall(PetscWeakFormGetResidual(ds->wf, NULL, 0, f, 0, &n0, &tmp0, &n1, &tmp1)); 1334 *f0 = tmp0 ? tmp0[0] : NULL; 1335 *f1 = tmp1 ? tmp1[0] : NULL; 1336 PetscFunctionReturn(PETSC_SUCCESS); 1337 } 1338 1339 /*@C 1340 PetscDSSetResidual - Set the pointwise residual function for a given test field 1341 1342 Not collective 1343 1344 Input Parameters: 1345 + ds - The `PetscDS` 1346 . f - The test field number 1347 . f0 - integrand for the test function term 1348 - f1 - integrand for the test function gradient term 1349 1350 Calling sequence for the callbacks f0 and f1: 1351 .vb 1352 f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1353 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1354 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1355 PetscReal t, const PetscReal x[], PetscScalar f0[]) 1356 .ve 1357 + dim - the spatial dimension 1358 . Nf - the number of fields 1359 . uOff - the offset into u[] and u_t[] for each field 1360 . uOff_x - the offset into u_x[] for each field 1361 . u - each field evaluated at the current point 1362 . u_t - the time derivative of each field evaluated at the current point 1363 . u_x - the gradient of each field evaluated at the current point 1364 . aOff - the offset into a[] and a_t[] for each auxiliary field 1365 . aOff_x - the offset into a_x[] for each auxiliary field 1366 . a - each auxiliary field evaluated at the current point 1367 . a_t - the time derivative of each auxiliary field evaluated at the current point 1368 . a_x - the gradient of auxiliary each field evaluated at the current point 1369 . t - current time 1370 . x - coordinates of the current point 1371 . numConstants - number of constant parameters 1372 . constants - constant parameters 1373 - f0 - output values at the current point 1374 1375 Level: intermediate 1376 1377 Note: 1378 We are using a first order FEM model for the weak form: \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t) 1379 1380 .seealso: `PetscDS`, `PetscDSGetResidual()` 1381 @*/ 1382 PetscErrorCode PetscDSSetResidual(PetscDS ds, PetscInt f, void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 1383 { 1384 PetscFunctionBegin; 1385 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1386 if (f0) PetscValidFunction(f0, 3); 1387 if (f1) PetscValidFunction(f1, 4); 1388 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 1389 PetscCall(PetscWeakFormSetIndexResidual(ds->wf, NULL, 0, f, 0, 0, f0, 0, f1)); 1390 PetscFunctionReturn(PETSC_SUCCESS); 1391 } 1392 1393 /*@C 1394 PetscDSGetRHSResidual - Get the pointwise RHS residual function for explicit timestepping for a given test field 1395 1396 Not collective 1397 1398 Input Parameters: 1399 + ds - The `PetscDS` 1400 - f - The test field number 1401 1402 Output Parameters: 1403 + f0 - integrand for the test function term 1404 - f1 - integrand for the test function gradient term 1405 1406 Calling sequence for the callbacks f0 and f1: 1407 .vb 1408 f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1409 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1410 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1411 PetscReal t, const PetscReal x[], PetscScalar f0[]) 1412 .ve 1413 + dim - the spatial dimension 1414 . Nf - the number of fields 1415 . uOff - the offset into u[] and u_t[] for each field 1416 . uOff_x - the offset into u_x[] for each field 1417 . u - each field evaluated at the current point 1418 . u_t - the time derivative of each field evaluated at the current point 1419 . u_x - the gradient of each field evaluated at the current point 1420 . aOff - the offset into a[] and a_t[] for each auxiliary field 1421 . aOff_x - the offset into a_x[] for each auxiliary field 1422 . a - each auxiliary field evaluated at the current point 1423 . a_t - the time derivative of each auxiliary field evaluated at the current point 1424 . a_x - the gradient of auxiliary each field evaluated at the current point 1425 . t - current time 1426 . x - coordinates of the current point 1427 . numConstants - number of constant parameters 1428 . constants - constant parameters 1429 - f0 - output values at the current point 1430 1431 Level: intermediate 1432 1433 Note: 1434 We are using a first order FEM model for the weak form: \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t) 1435 1436 .seealso: `PetscDS`, `PetscDSSetRHSResidual()` 1437 @*/ 1438 PetscErrorCode PetscDSGetRHSResidual(PetscDS ds, PetscInt f, void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 1439 { 1440 PetscPointFunc *tmp0, *tmp1; 1441 PetscInt n0, n1; 1442 1443 PetscFunctionBegin; 1444 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1445 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 1446 PetscCall(PetscWeakFormGetResidual(ds->wf, NULL, 0, f, 100, &n0, &tmp0, &n1, &tmp1)); 1447 *f0 = tmp0 ? tmp0[0] : NULL; 1448 *f1 = tmp1 ? tmp1[0] : NULL; 1449 PetscFunctionReturn(PETSC_SUCCESS); 1450 } 1451 1452 /*@C 1453 PetscDSSetRHSResidual - Set the pointwise residual function for explicit timestepping for a given test field 1454 1455 Not collective 1456 1457 Input Parameters: 1458 + ds - The `PetscDS` 1459 . f - The test field number 1460 . f0 - integrand for the test function term 1461 - f1 - integrand for the test function gradient term 1462 1463 Clling sequence for the callbacks f0 and f1: 1464 .vb 1465 f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1466 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1467 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1468 PetscReal t, const PetscReal x[], PetscScalar f0[]) 1469 .ve 1470 + dim - the spatial dimension 1471 . Nf - the number of fields 1472 . uOff - the offset into u[] and u_t[] for each field 1473 . uOff_x - the offset into u_x[] for each field 1474 . u - each field evaluated at the current point 1475 . u_t - the time derivative of each field evaluated at the current point 1476 . u_x - the gradient of each field evaluated at the current point 1477 . aOff - the offset into a[] and a_t[] for each auxiliary field 1478 . aOff_x - the offset into a_x[] for each auxiliary field 1479 . a - each auxiliary field evaluated at the current point 1480 . a_t - the time derivative of each auxiliary field evaluated at the current point 1481 . a_x - the gradient of auxiliary each field evaluated at the current point 1482 . t - current time 1483 . x - coordinates of the current point 1484 . numConstants - number of constant parameters 1485 . constants - constant parameters 1486 - f0 - output values at the current point 1487 1488 Level: intermediate 1489 1490 Note: 1491 We are using a first order FEM model for the weak form: \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t) 1492 1493 .seealso: `PetscDS`, `PetscDSGetResidual()` 1494 @*/ 1495 PetscErrorCode PetscDSSetRHSResidual(PetscDS ds, PetscInt f, void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 1496 { 1497 PetscFunctionBegin; 1498 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1499 if (f0) PetscValidFunction(f0, 3); 1500 if (f1) PetscValidFunction(f1, 4); 1501 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 1502 PetscCall(PetscWeakFormSetIndexResidual(ds->wf, NULL, 0, f, 100, 0, f0, 0, f1)); 1503 PetscFunctionReturn(PETSC_SUCCESS); 1504 } 1505 1506 /*@C 1507 PetscDSHasJacobian - Checks that the Jacobian functions have been set 1508 1509 Not collective 1510 1511 Input Parameter: 1512 . prob - The `PetscDS` 1513 1514 Output Parameter: 1515 . hasJac - flag that pointwise function for the Jacobian has been set 1516 1517 Level: intermediate 1518 1519 .seealso: `PetscDS`, `PetscDSGetJacobianPreconditioner()`, `PetscDSSetJacobianPreconditioner()`, `PetscDSGetJacobian()` 1520 @*/ 1521 PetscErrorCode PetscDSHasJacobian(PetscDS ds, PetscBool *hasJac) 1522 { 1523 PetscFunctionBegin; 1524 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1525 PetscCall(PetscWeakFormHasJacobian(ds->wf, hasJac)); 1526 PetscFunctionReturn(PETSC_SUCCESS); 1527 } 1528 1529 /*@C 1530 PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field 1531 1532 Not collective 1533 1534 Input Parameters: 1535 + ds - The `PetscDS` 1536 . f - The test field number 1537 - g - The field number 1538 1539 Output Parameters: 1540 + g0 - integrand for the test and basis function term 1541 . g1 - integrand for the test function and basis function gradient term 1542 . g2 - integrand for the test function gradient and basis function term 1543 - g3 - integrand for the test function gradient and basis function gradient term 1544 1545 Calling sequence for the callbacks g0, g1, g2 and g3: 1546 .vb 1547 g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1548 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1549 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1550 PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]) 1551 .ve 1552 + dim - the spatial dimension 1553 . Nf - the number of fields 1554 . uOff - the offset into u[] and u_t[] for each field 1555 . uOff_x - the offset into u_x[] for each field 1556 . u - each field evaluated at the current point 1557 . u_t - the time derivative of each field evaluated at the current point 1558 . u_x - the gradient of each field evaluated at the current point 1559 . aOff - the offset into a[] and a_t[] for each auxiliary field 1560 . aOff_x - the offset into a_x[] for each auxiliary field 1561 . a - each auxiliary field evaluated at the current point 1562 . a_t - the time derivative of each auxiliary field evaluated at the current point 1563 . a_x - the gradient of auxiliary each field evaluated at the current point 1564 . t - current time 1565 . u_tShift - the multiplier a for dF/dU_t 1566 . x - coordinates of the current point 1567 . numConstants - number of constant parameters 1568 . constants - constant parameters 1569 - g0 - output values at the current point 1570 1571 Level: intermediate 1572 1573 Note: 1574 We are using a first order FEM model for the weak form: 1575 \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi 1576 1577 .seealso: `PetscDS`, `PetscDSSetJacobian()` 1578 @*/ 1579 PetscErrorCode PetscDSGetJacobian(PetscDS ds, PetscInt f, PetscInt g, void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1580 { 1581 PetscPointJac *tmp0, *tmp1, *tmp2, *tmp3; 1582 PetscInt n0, n1, n2, n3; 1583 1584 PetscFunctionBegin; 1585 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1586 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 1587 PetscCheck(!(g < 0) && !(g >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", g, ds->Nf); 1588 PetscCall(PetscWeakFormGetJacobian(ds->wf, NULL, 0, f, g, 0, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3)); 1589 *g0 = tmp0 ? tmp0[0] : NULL; 1590 *g1 = tmp1 ? tmp1[0] : NULL; 1591 *g2 = tmp2 ? tmp2[0] : NULL; 1592 *g3 = tmp3 ? tmp3[0] : NULL; 1593 PetscFunctionReturn(PETSC_SUCCESS); 1594 } 1595 1596 /*@C 1597 PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields 1598 1599 Not collective 1600 1601 Input Parameters: 1602 + ds - The `PetscDS` 1603 . f - The test field number 1604 . g - The field number 1605 . g0 - integrand for the test and basis function term 1606 . g1 - integrand for the test function and basis function gradient term 1607 . g2 - integrand for the test function gradient and basis function term 1608 - g3 - integrand for the test function gradient and basis function gradient term 1609 1610 Calling sequence for the callbacks g0, g1, g2 and g3: 1611 .vb 1612 g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1613 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1614 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1615 PetscReal t, const PetscReal x[], PetscScalar g0[]) 1616 .ve 1617 + dim - the spatial dimension 1618 . Nf - the number of fields 1619 . uOff - the offset into u[] and u_t[] for each field 1620 . uOff_x - the offset into u_x[] for each field 1621 . u - each field evaluated at the current point 1622 . u_t - the time derivative of each field evaluated at the current point 1623 . u_x - the gradient of each field evaluated at the current point 1624 . aOff - the offset into a[] and a_t[] for each auxiliary field 1625 . aOff_x - the offset into a_x[] for each auxiliary field 1626 . a - each auxiliary field evaluated at the current point 1627 . a_t - the time derivative of each auxiliary field evaluated at the current point 1628 . a_x - the gradient of auxiliary each field evaluated at the current point 1629 . t - current time 1630 . u_tShift - the multiplier a for dF/dU_t 1631 . x - coordinates of the current point 1632 . numConstants - number of constant parameters 1633 . constants - constant parameters 1634 - g0 - output values at the current point 1635 1636 Level: intermediate 1637 1638 Note: 1639 We are using a first order FEM model for the weak form: 1640 \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi 1641 1642 .seealso: `PetscDS`, `PetscDSGetJacobian()` 1643 @*/ 1644 PetscErrorCode PetscDSSetJacobian(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1645 { 1646 PetscFunctionBegin; 1647 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1648 if (g0) PetscValidFunction(g0, 4); 1649 if (g1) PetscValidFunction(g1, 5); 1650 if (g2) PetscValidFunction(g2, 6); 1651 if (g3) PetscValidFunction(g3, 7); 1652 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 1653 PetscCheck(g >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", g); 1654 PetscCall(PetscWeakFormSetIndexJacobian(ds->wf, NULL, 0, f, g, 0, 0, g0, 0, g1, 0, g2, 0, g3)); 1655 PetscFunctionReturn(PETSC_SUCCESS); 1656 } 1657 1658 /*@C 1659 PetscDSUseJacobianPreconditioner - Set whether to construct a Jacobian preconditioner 1660 1661 Not collective 1662 1663 Input Parameters: 1664 + prob - The `PetscDS` 1665 - useJacPre - flag that enables construction of a Jacobian preconditioner 1666 1667 Level: intermediate 1668 1669 .seealso: `PetscDS`, `PetscDSGetJacobianPreconditioner()`, `PetscDSSetJacobianPreconditioner()`, `PetscDSGetJacobian()` 1670 @*/ 1671 PetscErrorCode PetscDSUseJacobianPreconditioner(PetscDS prob, PetscBool useJacPre) 1672 { 1673 PetscFunctionBegin; 1674 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1675 prob->useJacPre = useJacPre; 1676 PetscFunctionReturn(PETSC_SUCCESS); 1677 } 1678 1679 /*@C 1680 PetscDSHasJacobianPreconditioner - Checks if a Jacobian preconditioner matrix has been set 1681 1682 Not collective 1683 1684 Input Parameter: 1685 . prob - The `PetscDS` 1686 1687 Output Parameter: 1688 . hasJacPre - flag that pointwise function for Jacobian preconditioner matrix has been set 1689 1690 Level: intermediate 1691 1692 .seealso: `PetscDS`, `PetscDSGetJacobianPreconditioner()`, `PetscDSSetJacobianPreconditioner()`, `PetscDSGetJacobian()` 1693 @*/ 1694 PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS ds, PetscBool *hasJacPre) 1695 { 1696 PetscFunctionBegin; 1697 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1698 *hasJacPre = PETSC_FALSE; 1699 if (!ds->useJacPre) PetscFunctionReturn(PETSC_SUCCESS); 1700 PetscCall(PetscWeakFormHasJacobianPreconditioner(ds->wf, hasJacPre)); 1701 PetscFunctionReturn(PETSC_SUCCESS); 1702 } 1703 1704 /*@C 1705 PetscDSGetJacobianPreconditioner - Get the pointwise Jacobian preconditioner function for given test and basis field. If this is missing, 1706 the system matrix is used to build the preconditioner. 1707 1708 Not collective 1709 1710 Input Parameters: 1711 + ds - The `PetscDS` 1712 . f - The test field number 1713 - g - The field number 1714 1715 Output Parameters: 1716 + g0 - integrand for the test and basis function term 1717 . g1 - integrand for the test function and basis function gradient term 1718 . g2 - integrand for the test function gradient and basis function term 1719 - g3 - integrand for the test function gradient and basis function gradient term 1720 1721 Calling sequence for the callbacks g0, g1, g2 and g3: 1722 .vb 1723 g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1724 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1725 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1726 PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]) 1727 .ve 1728 + dim - the spatial dimension 1729 . Nf - the number of fields 1730 . uOff - the offset into u[] and u_t[] for each field 1731 . uOff_x - the offset into u_x[] for each field 1732 . u - each field evaluated at the current point 1733 . u_t - the time derivative of each field evaluated at the current point 1734 . u_x - the gradient of each field evaluated at the current point 1735 . aOff - the offset into a[] and a_t[] for each auxiliary field 1736 . aOff_x - the offset into a_x[] for each auxiliary field 1737 . a - each auxiliary field evaluated at the current point 1738 . a_t - the time derivative of each auxiliary field evaluated at the current point 1739 . a_x - the gradient of auxiliary each field evaluated at the current point 1740 . t - current time 1741 . u_tShift - the multiplier a for dF/dU_t 1742 . x - coordinates of the current point 1743 . numConstants - number of constant parameters 1744 . constants - constant parameters 1745 - g0 - output values at the current point 1746 1747 Level: intermediate 1748 1749 Note: 1750 We are using a first order FEM model for the weak form: 1751 \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi 1752 1753 .seealso: `PetscDS`, `PetscDSSetJacobianPreconditioner()`, `PetscDSGetJacobian()` 1754 @*/ 1755 PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1756 { 1757 PetscPointJac *tmp0, *tmp1, *tmp2, *tmp3; 1758 PetscInt n0, n1, n2, n3; 1759 1760 PetscFunctionBegin; 1761 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1762 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 1763 PetscCheck(!(g < 0) && !(g >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", g, ds->Nf); 1764 PetscCall(PetscWeakFormGetJacobianPreconditioner(ds->wf, NULL, 0, f, g, 0, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3)); 1765 *g0 = tmp0 ? tmp0[0] : NULL; 1766 *g1 = tmp1 ? tmp1[0] : NULL; 1767 *g2 = tmp2 ? tmp2[0] : NULL; 1768 *g3 = tmp3 ? tmp3[0] : NULL; 1769 PetscFunctionReturn(PETSC_SUCCESS); 1770 } 1771 1772 /*@C 1773 PetscDSSetJacobianPreconditioner - Set the pointwise Jacobian preconditioner function for given test and basis fields. 1774 If this is missing, the system matrix is used to build the preconditioner. 1775 1776 Not collective 1777 1778 Input Parameters: 1779 + ds - The `PetscDS` 1780 . f - The test field number 1781 . g - The field number 1782 . g0 - integrand for the test and basis function term 1783 . g1 - integrand for the test function and basis function gradient term 1784 . g2 - integrand for the test function gradient and basis function term 1785 - g3 - integrand for the test function gradient and basis function gradient term 1786 1787 Calling sequence for the callbacks g0, g1, g2 and g3: 1788 .vb 1789 g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1790 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1791 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1792 PetscReal t, const PetscReal x[], PetscScalar g0[]) 1793 .ve 1794 + dim - the spatial dimension 1795 . Nf - the number of fields 1796 . uOff - the offset into u[] and u_t[] for each field 1797 . uOff_x - the offset into u_x[] for each field 1798 . u - each field evaluated at the current point 1799 . u_t - the time derivative of each field evaluated at the current point 1800 . u_x - the gradient of each field evaluated at the current point 1801 . aOff - the offset into a[] and a_t[] for each auxiliary field 1802 . aOff_x - the offset into a_x[] for each auxiliary field 1803 . a - each auxiliary field evaluated at the current point 1804 . a_t - the time derivative of each auxiliary field evaluated at the current point 1805 . a_x - the gradient of auxiliary each field evaluated at the current point 1806 . t - current time 1807 . u_tShift - the multiplier a for dF/dU_t 1808 . x - coordinates of the current point 1809 . numConstants - number of constant parameters 1810 . constants - constant parameters 1811 - g0 - output values at the current point 1812 1813 Level: intermediate 1814 1815 Note: 1816 We are using a first order FEM model for the weak form: 1817 \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi 1818 1819 .seealso: `PetscDS`, `PetscDSGetJacobianPreconditioner()`, `PetscDSSetJacobian()` 1820 @*/ 1821 PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1822 { 1823 PetscFunctionBegin; 1824 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1825 if (g0) PetscValidFunction(g0, 4); 1826 if (g1) PetscValidFunction(g1, 5); 1827 if (g2) PetscValidFunction(g2, 6); 1828 if (g3) PetscValidFunction(g3, 7); 1829 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 1830 PetscCheck(g >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", g); 1831 PetscCall(PetscWeakFormSetIndexJacobianPreconditioner(ds->wf, NULL, 0, f, g, 0, 0, g0, 0, g1, 0, g2, 0, g3)); 1832 PetscFunctionReturn(PETSC_SUCCESS); 1833 } 1834 1835 /*@C 1836 PetscDSHasDynamicJacobian - Signals that a dynamic Jacobian, dF/du_t, has been set 1837 1838 Not collective 1839 1840 Input Parameter: 1841 . ds - The `PetscDS` 1842 1843 Output Parameter: 1844 . hasDynJac - flag that pointwise function for dynamic Jacobian has been set 1845 1846 Level: intermediate 1847 1848 .seealso: `PetscDS`, `PetscDSGetDynamicJacobian()`, `PetscDSSetDynamicJacobian()`, `PetscDSGetJacobian()` 1849 @*/ 1850 PetscErrorCode PetscDSHasDynamicJacobian(PetscDS ds, PetscBool *hasDynJac) 1851 { 1852 PetscFunctionBegin; 1853 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1854 PetscCall(PetscWeakFormHasDynamicJacobian(ds->wf, hasDynJac)); 1855 PetscFunctionReturn(PETSC_SUCCESS); 1856 } 1857 1858 /*@C 1859 PetscDSGetDynamicJacobian - Get the pointwise dynamic Jacobian, dF/du_t, function for given test and basis field 1860 1861 Not collective 1862 1863 Input Parameters: 1864 + ds - The `PetscDS` 1865 . f - The test field number 1866 - g - The field number 1867 1868 Output Parameters: 1869 + g0 - integrand for the test and basis function term 1870 . g1 - integrand for the test function and basis function gradient term 1871 . g2 - integrand for the test function gradient and basis function term 1872 - g3 - integrand for the test function gradient and basis function gradient term 1873 1874 Calling sequence for the callbacks g0, g1, g2 and g3: 1875 .vb 1876 g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1877 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1878 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1879 PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]) 1880 .ve 1881 + dim - the spatial dimension 1882 . Nf - the number of fields 1883 . uOff - the offset into u[] and u_t[] for each field 1884 . uOff_x - the offset into u_x[] for each field 1885 . u - each field evaluated at the current point 1886 . u_t - the time derivative of each field evaluated at the current point 1887 . u_x - the gradient of each field evaluated at the current point 1888 . aOff - the offset into a[] and a_t[] for each auxiliary field 1889 . aOff_x - the offset into a_x[] for each auxiliary field 1890 . a - each auxiliary field evaluated at the current point 1891 . a_t - the time derivative of each auxiliary field evaluated at the current point 1892 . a_x - the gradient of auxiliary each field evaluated at the current point 1893 . t - current time 1894 . u_tShift - the multiplier a for dF/dU_t 1895 . x - coordinates of the current point 1896 . numConstants - number of constant parameters 1897 . constants - constant parameters 1898 - g0 - output values at the current point 1899 1900 Level: intermediate 1901 1902 Note: 1903 We are using a first order FEM model for the weak form: 1904 \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi 1905 1906 .seealso: `PetscDS`, `PetscDSSetJacobian()` 1907 @*/ 1908 PetscErrorCode PetscDSGetDynamicJacobian(PetscDS ds, PetscInt f, PetscInt g, void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1909 { 1910 PetscPointJac *tmp0, *tmp1, *tmp2, *tmp3; 1911 PetscInt n0, n1, n2, n3; 1912 1913 PetscFunctionBegin; 1914 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1915 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 1916 PetscCheck(!(g < 0) && !(g >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", g, ds->Nf); 1917 PetscCall(PetscWeakFormGetDynamicJacobian(ds->wf, NULL, 0, f, g, 0, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3)); 1918 *g0 = tmp0 ? tmp0[0] : NULL; 1919 *g1 = tmp1 ? tmp1[0] : NULL; 1920 *g2 = tmp2 ? tmp2[0] : NULL; 1921 *g3 = tmp3 ? tmp3[0] : NULL; 1922 PetscFunctionReturn(PETSC_SUCCESS); 1923 } 1924 1925 /*@C 1926 PetscDSSetDynamicJacobian - Set the pointwise dynamic Jacobian, dF/du_t, function for given test and basis fields 1927 1928 Not collective 1929 1930 Input Parameters: 1931 + ds - The `PetscDS` 1932 . f - The test field number 1933 . g - The field number 1934 . g0 - integrand for the test and basis function term 1935 . g1 - integrand for the test function and basis function gradient term 1936 . g2 - integrand for the test function gradient and basis function term 1937 - g3 - integrand for the test function gradient and basis function gradient term 1938 1939 Calling sequence for the callbacks g0, g1, g2 and g3: 1940 .vb 1941 g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1942 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1943 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1944 PetscReal t, const PetscReal x[], PetscScalar g0[]) 1945 .ve 1946 + dim - the spatial dimension 1947 . Nf - the number of fields 1948 . uOff - the offset into u[] and u_t[] for each field 1949 . uOff_x - the offset into u_x[] for each field 1950 . u - each field evaluated at the current point 1951 . u_t - the time derivative of each field evaluated at the current point 1952 . u_x - the gradient of each field evaluated at the current point 1953 . aOff - the offset into a[] and a_t[] for each auxiliary field 1954 . aOff_x - the offset into a_x[] for each auxiliary field 1955 . a - each auxiliary field evaluated at the current point 1956 . a_t - the time derivative of each auxiliary field evaluated at the current point 1957 . a_x - the gradient of auxiliary each field evaluated at the current point 1958 . t - current time 1959 . u_tShift - the multiplier a for dF/dU_t 1960 . x - coordinates of the current point 1961 . numConstants - number of constant parameters 1962 . constants - constant parameters 1963 - g0 - output values at the current point 1964 1965 Level: intermediate 1966 1967 Note: 1968 We are using a first order FEM model for the weak form: 1969 \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi 1970 1971 .seealso: `PetscDS`, `PetscDSGetJacobian()` 1972 @*/ 1973 PetscErrorCode PetscDSSetDynamicJacobian(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1974 { 1975 PetscFunctionBegin; 1976 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1977 if (g0) PetscValidFunction(g0, 4); 1978 if (g1) PetscValidFunction(g1, 5); 1979 if (g2) PetscValidFunction(g2, 6); 1980 if (g3) PetscValidFunction(g3, 7); 1981 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 1982 PetscCheck(g >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", g); 1983 PetscCall(PetscWeakFormSetIndexDynamicJacobian(ds->wf, NULL, 0, f, g, 0, 0, g0, 0, g1, 0, g2, 0, g3)); 1984 PetscFunctionReturn(PETSC_SUCCESS); 1985 } 1986 1987 /*@C 1988 PetscDSGetRiemannSolver - Returns the Riemann solver for the given field 1989 1990 Not collective 1991 1992 Input Parameters: 1993 + ds - The `PetscDS` object 1994 - f - The field number 1995 1996 Output Parameter: 1997 . r - Riemann solver 1998 1999 Calling sequence for r: 2000 .vb 2001 r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 2002 .ve 2003 + dim - The spatial dimension 2004 . Nf - The number of fields 2005 . x - The coordinates at a point on the interface 2006 . n - The normal vector to the interface 2007 . uL - The state vector to the left of the interface 2008 . uR - The state vector to the right of the interface 2009 . flux - output array of flux through the interface 2010 . numConstants - number of constant parameters 2011 . constants - constant parameters 2012 - ctx - optional user context 2013 2014 Level: intermediate 2015 2016 .seealso: `PetscDS`, `PetscDSSetRiemannSolver()` 2017 @*/ 2018 PetscErrorCode PetscDSGetRiemannSolver(PetscDS ds, PetscInt f, void (**r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscInt numConstants, const PetscScalar constants[], PetscScalar flux[], void *ctx)) 2019 { 2020 PetscRiemannFunc *tmp; 2021 PetscInt n; 2022 2023 PetscFunctionBegin; 2024 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2025 PetscValidPointer(r, 3); 2026 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 2027 PetscCall(PetscWeakFormGetRiemannSolver(ds->wf, NULL, 0, f, 0, &n, &tmp)); 2028 *r = tmp ? tmp[0] : NULL; 2029 PetscFunctionReturn(PETSC_SUCCESS); 2030 } 2031 2032 /*@C 2033 PetscDSSetRiemannSolver - Sets the Riemann solver for the given field 2034 2035 Not collective 2036 2037 Input Parameters: 2038 + ds - The `PetscDS` object 2039 . f - The field number 2040 - r - Riemann solver 2041 2042 Calling sequence for r: 2043 .vb 2044 r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 2045 .ve 2046 + dim - The spatial dimension 2047 . Nf - The number of fields 2048 . x - The coordinates at a point on the interface 2049 . n - The normal vector to the interface 2050 . uL - The state vector to the left of the interface 2051 . uR - The state vector to the right of the interface 2052 . flux - output array of flux through the interface 2053 . numConstants - number of constant parameters 2054 . constants - constant parameters 2055 - ctx - optional user context 2056 2057 Level: intermediate 2058 2059 .seealso: `PetscDS`, `PetscDSGetRiemannSolver()` 2060 @*/ 2061 PetscErrorCode PetscDSSetRiemannSolver(PetscDS ds, PetscInt f, void (*r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscInt numConstants, const PetscScalar constants[], PetscScalar flux[], void *ctx)) 2062 { 2063 PetscFunctionBegin; 2064 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2065 if (r) PetscValidFunction(r, 3); 2066 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 2067 PetscCall(PetscWeakFormSetIndexRiemannSolver(ds->wf, NULL, 0, f, 0, 0, r)); 2068 PetscFunctionReturn(PETSC_SUCCESS); 2069 } 2070 2071 /*@C 2072 PetscDSGetUpdate - Get the pointwise update function for a given field 2073 2074 Not collective 2075 2076 Input Parameters: 2077 + ds - The `PetscDS` 2078 - f - The field number 2079 2080 Output Parameter: 2081 . update - update function 2082 2083 Calling sequence for the callback update: 2084 .vb 2085 update(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2086 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2087 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2088 PetscReal t, const PetscReal x[], PetscScalar uNew[]) 2089 .ve 2090 + dim - the spatial dimension 2091 . Nf - the number of fields 2092 . uOff - the offset into u[] and u_t[] for each field 2093 . uOff_x - the offset into u_x[] for each field 2094 . u - each field evaluated at the current point 2095 . u_t - the time derivative of each field evaluated at the current point 2096 . u_x - the gradient of each field evaluated at the current point 2097 . aOff - the offset into a[] and a_t[] for each auxiliary field 2098 . aOff_x - the offset into a_x[] for each auxiliary field 2099 . a - each auxiliary field evaluated at the current point 2100 . a_t - the time derivative of each auxiliary field evaluated at the current point 2101 . a_x - the gradient of auxiliary each field evaluated at the current point 2102 . t - current time 2103 . x - coordinates of the current point 2104 - uNew - new value for field at the current point 2105 2106 Level: intermediate 2107 2108 .seealso: `PetscDS`, `PetscDSSetUpdate()`, `PetscDSSetResidual()` 2109 @*/ 2110 PetscErrorCode PetscDSGetUpdate(PetscDS ds, PetscInt f, void (**update)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uNew[])) 2111 { 2112 PetscFunctionBegin; 2113 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2114 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 2115 if (update) { 2116 PetscValidPointer(update, 3); 2117 *update = ds->update[f]; 2118 } 2119 PetscFunctionReturn(PETSC_SUCCESS); 2120 } 2121 2122 /*@C 2123 PetscDSSetUpdate - Set the pointwise update function for a given field 2124 2125 Not collective 2126 2127 Input Parameters: 2128 + ds - The `PetscDS` 2129 . f - The field number 2130 - update - update function 2131 2132 Calling sequence for the callback update: 2133 .vb 2134 update(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2135 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2136 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2137 PetscReal t, const PetscReal x[], PetscScalar uNew[]) 2138 .ve 2139 + dim - the spatial dimension 2140 . Nf - the number of fields 2141 . uOff - the offset into u[] and u_t[] for each field 2142 . uOff_x - the offset into u_x[] for each field 2143 . u - each field evaluated at the current point 2144 . u_t - the time derivative of each field evaluated at the current point 2145 . u_x - the gradient of each field evaluated at the current point 2146 . aOff - the offset into a[] and a_t[] for each auxiliary field 2147 . aOff_x - the offset into a_x[] for each auxiliary field 2148 . a - each auxiliary field evaluated at the current point 2149 . a_t - the time derivative of each auxiliary field evaluated at the current point 2150 . a_x - the gradient of auxiliary each field evaluated at the current point 2151 . t - current time 2152 . x - coordinates of the current point 2153 - uNew - new field values at the current point 2154 2155 Level: intermediate 2156 2157 .seealso: `PetscDS`, `PetscDSGetResidual()` 2158 @*/ 2159 PetscErrorCode PetscDSSetUpdate(PetscDS ds, PetscInt f, void (*update)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uNew[])) 2160 { 2161 PetscFunctionBegin; 2162 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2163 if (update) PetscValidFunction(update, 3); 2164 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 2165 PetscCall(PetscDSEnlarge_Static(ds, f + 1)); 2166 ds->update[f] = update; 2167 PetscFunctionReturn(PETSC_SUCCESS); 2168 } 2169 2170 PetscErrorCode PetscDSGetContext(PetscDS ds, PetscInt f, void *ctx) 2171 { 2172 PetscFunctionBegin; 2173 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2174 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 2175 PetscValidPointer(ctx, 3); 2176 *(void **)ctx = ds->ctx[f]; 2177 PetscFunctionReturn(PETSC_SUCCESS); 2178 } 2179 2180 PetscErrorCode PetscDSSetContext(PetscDS ds, PetscInt f, void *ctx) 2181 { 2182 PetscFunctionBegin; 2183 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2184 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 2185 PetscCall(PetscDSEnlarge_Static(ds, f + 1)); 2186 ds->ctx[f] = ctx; 2187 PetscFunctionReturn(PETSC_SUCCESS); 2188 } 2189 2190 /*@C 2191 PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field 2192 2193 Not collective 2194 2195 Input Parameters: 2196 + ds - The PetscDS 2197 - f - The test field number 2198 2199 Output Parameters: 2200 + f0 - boundary integrand for the test function term 2201 - f1 - boundary integrand for the test function gradient term 2202 2203 Calling sequence for the callbacks f0 and f1: 2204 .vb 2205 f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2206 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2207 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2208 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]) 2209 .ve 2210 + dim - the spatial dimension 2211 . Nf - the number of fields 2212 . uOff - the offset into u[] and u_t[] for each field 2213 . uOff_x - the offset into u_x[] for each field 2214 . u - each field evaluated at the current point 2215 . u_t - the time derivative of each field evaluated at the current point 2216 . u_x - the gradient of each field evaluated at the current point 2217 . aOff - the offset into a[] and a_t[] for each auxiliary field 2218 . aOff_x - the offset into a_x[] for each auxiliary field 2219 . a - each auxiliary field evaluated at the current point 2220 . a_t - the time derivative of each auxiliary field evaluated at the current point 2221 . a_x - the gradient of auxiliary each field evaluated at the current point 2222 . t - current time 2223 . x - coordinates of the current point 2224 . n - unit normal at the current point 2225 . numConstants - number of constant parameters 2226 . constants - constant parameters 2227 - f0 - output values at the current point 2228 2229 Level: intermediate 2230 2231 Note: 2232 We are using a first order FEM model for the weak form: 2233 \int_\Gamma \phi {\vec f}_0(u, u_t, \nabla u, x, t) \cdot \hat n + \nabla\phi \cdot {\overleftrightarrow f}_1(u, u_t, \nabla u, x, t) \cdot \hat n 2234 2235 .seealso: `PetscDS`, `PetscDSSetBdResidual()` 2236 @*/ 2237 PetscErrorCode PetscDSGetBdResidual(PetscDS ds, PetscInt f, void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 2238 { 2239 PetscBdPointFunc *tmp0, *tmp1; 2240 PetscInt n0, n1; 2241 2242 PetscFunctionBegin; 2243 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2244 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 2245 PetscCall(PetscWeakFormGetBdResidual(ds->wf, NULL, 0, f, 0, &n0, &tmp0, &n1, &tmp1)); 2246 *f0 = tmp0 ? tmp0[0] : NULL; 2247 *f1 = tmp1 ? tmp1[0] : NULL; 2248 PetscFunctionReturn(PETSC_SUCCESS); 2249 } 2250 2251 /*@C 2252 PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field 2253 2254 Not collective 2255 2256 Input Parameters: 2257 + ds - The `PetscDS` 2258 . f - The test field number 2259 . f0 - boundary integrand for the test function term 2260 - f1 - boundary integrand for the test function gradient term 2261 2262 Calling sequence for the callbacks f0 and f1: 2263 .vb 2264 f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2265 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2266 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2267 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]) 2268 .ve 2269 + dim - the spatial dimension 2270 . Nf - the number of fields 2271 . uOff - the offset into u[] and u_t[] for each field 2272 . uOff_x - the offset into u_x[] for each field 2273 . u - each field evaluated at the current point 2274 . u_t - the time derivative of each field evaluated at the current point 2275 . u_x - the gradient of each field evaluated at the current point 2276 . aOff - the offset into a[] and a_t[] for each auxiliary field 2277 . aOff_x - the offset into a_x[] for each auxiliary field 2278 . a - each auxiliary field evaluated at the current point 2279 . a_t - the time derivative of each auxiliary field evaluated at the current point 2280 . a_x - the gradient of auxiliary each field evaluated at the current point 2281 . t - current time 2282 . x - coordinates of the current point 2283 . n - unit normal at the current point 2284 . numConstants - number of constant parameters 2285 . constants - constant parameters 2286 - f0 - output values at the current point 2287 2288 Level: intermediate 2289 2290 Note: 2291 We are using a first order FEM model for the weak form: 2292 \int_\Gamma \phi {\vec f}_0(u, u_t, \nabla u, x, t) \cdot \hat n + \nabla\phi \cdot {\overleftrightarrow f}_1(u, u_t, \nabla u, x, t) \cdot \hat n 2293 2294 .seealso: `PetscDS`, `PetscDSGetBdResidual()` 2295 @*/ 2296 PetscErrorCode PetscDSSetBdResidual(PetscDS ds, PetscInt f, void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 2297 { 2298 PetscFunctionBegin; 2299 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2300 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 2301 PetscCall(PetscWeakFormSetIndexBdResidual(ds->wf, NULL, 0, f, 0, 0, f0, 0, f1)); 2302 PetscFunctionReturn(PETSC_SUCCESS); 2303 } 2304 2305 /*@ 2306 PetscDSHasBdJacobian - Indicates that boundary Jacobian functions have been set 2307 2308 Not collective 2309 2310 Input Parameter: 2311 . ds - The `PetscDS` 2312 2313 Output Parameter: 2314 . hasBdJac - flag that pointwise function for the boundary Jacobian has been set 2315 2316 Level: intermediate 2317 2318 .seealso: `PetscDS`, `PetscDSHasJacobian()`, `PetscDSSetBdJacobian()`, `PetscDSGetBdJacobian()` 2319 @*/ 2320 PetscErrorCode PetscDSHasBdJacobian(PetscDS ds, PetscBool *hasBdJac) 2321 { 2322 PetscFunctionBegin; 2323 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2324 PetscValidBoolPointer(hasBdJac, 2); 2325 PetscCall(PetscWeakFormHasBdJacobian(ds->wf, hasBdJac)); 2326 PetscFunctionReturn(PETSC_SUCCESS); 2327 } 2328 2329 /*@C 2330 PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field 2331 2332 Not collective 2333 2334 Input Parameters: 2335 + ds - The `PetscDS` 2336 . f - The test field number 2337 - g - The field number 2338 2339 Output Parameters: 2340 + g0 - integrand for the test and basis function term 2341 . g1 - integrand for the test function and basis function gradient term 2342 . g2 - integrand for the test function gradient and basis function term 2343 - g3 - integrand for the test function gradient and basis function gradient term 2344 2345 Calling sequence for the callbacks g0, g1, g2 and g3: 2346 .vb 2347 g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2348 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2349 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2350 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[]) 2351 .ve 2352 + dim - the spatial dimension 2353 . Nf - the number of fields 2354 . uOff - the offset into u[] and u_t[] for each field 2355 . uOff_x - the offset into u_x[] for each field 2356 . u - each field evaluated at the current point 2357 . u_t - the time derivative of each field evaluated at the current point 2358 . u_x - the gradient of each field evaluated at the current point 2359 . aOff - the offset into a[] and a_t[] for each auxiliary field 2360 . aOff_x - the offset into a_x[] for each auxiliary field 2361 . a - each auxiliary field evaluated at the current point 2362 . a_t - the time derivative of each auxiliary field evaluated at the current point 2363 . a_x - the gradient of auxiliary each field evaluated at the current point 2364 . t - current time 2365 . u_tShift - the multiplier a for dF/dU_t 2366 . x - coordinates of the current point 2367 . n - normal at the current point 2368 . numConstants - number of constant parameters 2369 . constants - constant parameters 2370 - g0 - output values at the current point 2371 2372 Level: intermediate 2373 2374 Note: 2375 We are using a first order FEM model for the weak form: 2376 \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi 2377 2378 .seealso: `PetscDS`, `PetscDSSetBdJacobian()` 2379 @*/ 2380 PetscErrorCode PetscDSGetBdJacobian(PetscDS ds, PetscInt f, PetscInt g, void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 2381 { 2382 PetscBdPointJac *tmp0, *tmp1, *tmp2, *tmp3; 2383 PetscInt n0, n1, n2, n3; 2384 2385 PetscFunctionBegin; 2386 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2387 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 2388 PetscCheck(!(g < 0) && !(g >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", g, ds->Nf); 2389 PetscCall(PetscWeakFormGetBdJacobian(ds->wf, NULL, 0, f, g, 0, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3)); 2390 *g0 = tmp0 ? tmp0[0] : NULL; 2391 *g1 = tmp1 ? tmp1[0] : NULL; 2392 *g2 = tmp2 ? tmp2[0] : NULL; 2393 *g3 = tmp3 ? tmp3[0] : NULL; 2394 PetscFunctionReturn(PETSC_SUCCESS); 2395 } 2396 2397 /*@C 2398 PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field 2399 2400 Not collective 2401 2402 Input Parameters: 2403 + ds - The PetscDS 2404 . f - The test field number 2405 . g - The field number 2406 . g0 - integrand for the test and basis function term 2407 . g1 - integrand for the test function and basis function gradient term 2408 . g2 - integrand for the test function gradient and basis function term 2409 - g3 - integrand for the test function gradient and basis function gradient term 2410 2411 Calling sequence for the callbacks g0, g1, g2 and g3: 2412 .vb 2413 g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2414 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2415 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2416 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[]) 2417 .ve 2418 + dim - the spatial dimension 2419 . Nf - the number of fields 2420 . uOff - the offset into u[] and u_t[] for each field 2421 . uOff_x - the offset into u_x[] for each field 2422 . u - each field evaluated at the current point 2423 . u_t - the time derivative of each field evaluated at the current point 2424 . u_x - the gradient of each field evaluated at the current point 2425 . aOff - the offset into a[] and a_t[] for each auxiliary field 2426 . aOff_x - the offset into a_x[] for each auxiliary field 2427 . a - each auxiliary field evaluated at the current point 2428 . a_t - the time derivative of each auxiliary field evaluated at the current point 2429 . a_x - the gradient of auxiliary each field evaluated at the current point 2430 . t - current time 2431 . u_tShift - the multiplier a for dF/dU_t 2432 . x - coordinates of the current point 2433 . n - normal at the current point 2434 . numConstants - number of constant parameters 2435 . constants - constant parameters 2436 - g0 - output values at the current point 2437 2438 Level: intermediate 2439 2440 Note: 2441 We are using a first order FEM model for the weak form: 2442 \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi 2443 2444 .seealso: `PetscDS`, `PetscDSGetBdJacobian()` 2445 @*/ 2446 PetscErrorCode PetscDSSetBdJacobian(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 2447 { 2448 PetscFunctionBegin; 2449 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2450 if (g0) PetscValidFunction(g0, 4); 2451 if (g1) PetscValidFunction(g1, 5); 2452 if (g2) PetscValidFunction(g2, 6); 2453 if (g3) PetscValidFunction(g3, 7); 2454 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 2455 PetscCheck(g >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", g); 2456 PetscCall(PetscWeakFormSetIndexBdJacobian(ds->wf, NULL, 0, f, g, 0, 0, g0, 0, g1, 0, g2, 0, g3)); 2457 PetscFunctionReturn(PETSC_SUCCESS); 2458 } 2459 2460 /*@ 2461 PetscDSHasBdJacobianPreconditioner - Signals that boundary Jacobian preconditioner functions have been set 2462 2463 Not collective 2464 2465 Input Parameter: 2466 . ds - The `PetscDS` 2467 2468 Output Parameter: 2469 . hasBdJac - flag that pointwise function for the boundary Jacobian preconditioner has been set 2470 2471 Level: intermediate 2472 2473 .seealso: `PetscDS`, `PetscDSHasJacobian()`, `PetscDSSetBdJacobian()`, `PetscDSGetBdJacobian()` 2474 @*/ 2475 PetscErrorCode PetscDSHasBdJacobianPreconditioner(PetscDS ds, PetscBool *hasBdJacPre) 2476 { 2477 PetscFunctionBegin; 2478 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2479 PetscValidBoolPointer(hasBdJacPre, 2); 2480 PetscCall(PetscWeakFormHasBdJacobianPreconditioner(ds->wf, hasBdJacPre)); 2481 PetscFunctionReturn(PETSC_SUCCESS); 2482 } 2483 2484 /*@C 2485 PetscDSGetBdJacobianPreconditioner - Get the pointwise boundary Jacobian preconditioner function for given test and basis field 2486 2487 Not collective 2488 2489 Input Parameters: 2490 + ds - The `PetscDS` 2491 . f - The test field number 2492 - g - The field number 2493 2494 Output Parameters: 2495 + g0 - integrand for the test and basis function term 2496 . g1 - integrand for the test function and basis function gradient term 2497 . g2 - integrand for the test function gradient and basis function term 2498 - g3 - integrand for the test function gradient and basis function gradient term 2499 2500 Calling sequence for the callbacks g0, g1, g2 and g3: 2501 .vb 2502 g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2503 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2504 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2505 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 2506 .ve 2507 + dim - the spatial dimension 2508 . Nf - the number of fields 2509 . NfAux - the number of auxiliary fields 2510 . uOff - the offset into u[] and u_t[] for each field 2511 . uOff_x - the offset into u_x[] for each field 2512 . u - each field evaluated at the current point 2513 . u_t - the time derivative of each field evaluated at the current point 2514 . u_x - the gradient of each field evaluated at the current point 2515 . aOff - the offset into a[] and a_t[] for each auxiliary field 2516 . aOff_x - the offset into a_x[] for each auxiliary field 2517 . a - each auxiliary field evaluated at the current point 2518 . a_t - the time derivative of each auxiliary field evaluated at the current point 2519 . a_x - the gradient of auxiliary each field evaluated at the current point 2520 . t - current time 2521 . u_tShift - the multiplier a for dF/dU_t 2522 . x - coordinates of the current point 2523 . n - normal at the current point 2524 . numConstants - number of constant parameters 2525 . constants - constant parameters 2526 - g0 - output values at the current point 2527 2528 Level: intermediate 2529 2530 Note: 2531 We are using a first order FEM model for the weak form: 2532 \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi 2533 2534 Fortran Note: 2535 This is not yet available in Fortran. 2536 2537 .seealso: `PetscDS`, `PetscDSSetBdJacobianPreconditioner()` 2538 @*/ 2539 PetscErrorCode PetscDSGetBdJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 2540 { 2541 PetscBdPointJac *tmp0, *tmp1, *tmp2, *tmp3; 2542 PetscInt n0, n1, n2, n3; 2543 2544 PetscFunctionBegin; 2545 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2546 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 2547 PetscCheck(!(g < 0) && !(g >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", g, ds->Nf); 2548 PetscCall(PetscWeakFormGetBdJacobianPreconditioner(ds->wf, NULL, 0, f, g, 0, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3)); 2549 *g0 = tmp0 ? tmp0[0] : NULL; 2550 *g1 = tmp1 ? tmp1[0] : NULL; 2551 *g2 = tmp2 ? tmp2[0] : NULL; 2552 *g3 = tmp3 ? tmp3[0] : NULL; 2553 PetscFunctionReturn(PETSC_SUCCESS); 2554 } 2555 2556 /*@C 2557 PetscDSSetBdJacobianPreconditioner - Set the pointwise boundary Jacobian preconditioner function for given test and basis field 2558 2559 Not collective 2560 2561 Input Parameters: 2562 + ds - The `PetscDS` 2563 . f - The test field number 2564 . g - The field number 2565 . g0 - integrand for the test and basis function term 2566 . g1 - integrand for the test function and basis function gradient term 2567 . g2 - integrand for the test function gradient and basis function term 2568 - g3 - integrand for the test function gradient and basis function gradient term 2569 2570 Calling sequence for the callbacks g0, g1, g2 and g3: 2571 .vb 2572 g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2573 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2574 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2575 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 2576 .ve 2577 + dim - the spatial dimension 2578 . Nf - the number of fields 2579 . NfAux - the number of auxiliary fields 2580 . uOff - the offset into u[] and u_t[] for each field 2581 . uOff_x - the offset into u_x[] for each field 2582 . u - each field evaluated at the current point 2583 . u_t - the time derivative of each field evaluated at the current point 2584 . u_x - the gradient of each field evaluated at the current point 2585 . aOff - the offset into a[] and a_t[] for each auxiliary field 2586 . aOff_x - the offset into a_x[] for each auxiliary field 2587 . a - each auxiliary field evaluated at the current point 2588 . a_t - the time derivative of each auxiliary field evaluated at the current point 2589 . a_x - the gradient of auxiliary each field evaluated at the current point 2590 . t - current time 2591 . u_tShift - the multiplier a for dF/dU_t 2592 . x - coordinates of the current point 2593 . n - normal at the current point 2594 . numConstants - number of constant parameters 2595 . constants - constant parameters 2596 - g0 - output values at the current point 2597 2598 Level: intermediate 2599 2600 Note: 2601 We are using a first order FEM model for the weak form: 2602 \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi 2603 2604 Fortran Note: 2605 This is not yet available in Fortran. 2606 2607 .seealso: `PetscDS`, `PetscDSGetBdJacobianPreconditioner()` 2608 @*/ 2609 PetscErrorCode PetscDSSetBdJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 2610 { 2611 PetscFunctionBegin; 2612 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2613 if (g0) PetscValidFunction(g0, 4); 2614 if (g1) PetscValidFunction(g1, 5); 2615 if (g2) PetscValidFunction(g2, 6); 2616 if (g3) PetscValidFunction(g3, 7); 2617 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 2618 PetscCheck(g >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", g); 2619 PetscCall(PetscWeakFormSetIndexBdJacobianPreconditioner(ds->wf, NULL, 0, f, g, 0, 0, g0, 0, g1, 0, g2, 0, g3)); 2620 PetscFunctionReturn(PETSC_SUCCESS); 2621 } 2622 2623 /*@C 2624 PetscDSGetExactSolution - Get the pointwise exact solution function for a given test field 2625 2626 Not collective 2627 2628 Input Parameters: 2629 + prob - The PetscDS 2630 - f - The test field number 2631 2632 Output Parameters: 2633 + exactSol - exact solution for the test field 2634 - exactCtx - exact solution context 2635 2636 Calling sequence for the solution functions: 2637 .vb 2638 sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) 2639 .ve 2640 + dim - the spatial dimension 2641 . t - current time 2642 . x - coordinates of the current point 2643 . Nc - the number of field components 2644 . u - the solution field evaluated at the current point 2645 - ctx - a user context 2646 2647 Level: intermediate 2648 2649 .seealso: `PetscDS`, `PetscDSSetExactSolution()`, `PetscDSGetExactSolutionTimeDerivative()` 2650 @*/ 2651 PetscErrorCode PetscDSGetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (**sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void **ctx) 2652 { 2653 PetscFunctionBegin; 2654 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2655 PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf); 2656 if (sol) { 2657 PetscValidPointer(sol, 3); 2658 *sol = prob->exactSol[f]; 2659 } 2660 if (ctx) { 2661 PetscValidPointer(ctx, 4); 2662 *ctx = prob->exactCtx[f]; 2663 } 2664 PetscFunctionReturn(PETSC_SUCCESS); 2665 } 2666 2667 /*@C 2668 PetscDSSetExactSolution - Set the pointwise exact solution function for a given test field 2669 2670 Not collective 2671 2672 Input Parameters: 2673 + prob - The `PetscDS` 2674 . f - The test field number 2675 . sol - solution function for the test fields 2676 - ctx - solution context or NULL 2677 2678 Calling sequence for solution functions: 2679 .vb 2680 sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) 2681 .ve 2682 + dim - the spatial dimension 2683 . t - current time 2684 . x - coordinates of the current point 2685 . Nc - the number of field components 2686 . u - the solution field evaluated at the current point 2687 - ctx - a user context 2688 2689 Level: intermediate 2690 2691 .seealso: `PetscDS`, `PetscDSGetExactSolution()` 2692 @*/ 2693 PetscErrorCode PetscDSSetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (*sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void *ctx) 2694 { 2695 PetscFunctionBegin; 2696 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2697 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 2698 PetscCall(PetscDSEnlarge_Static(prob, f + 1)); 2699 if (sol) { 2700 PetscValidFunction(sol, 3); 2701 prob->exactSol[f] = sol; 2702 } 2703 if (ctx) { 2704 PetscValidFunction(ctx, 4); 2705 prob->exactCtx[f] = ctx; 2706 } 2707 PetscFunctionReturn(PETSC_SUCCESS); 2708 } 2709 2710 /*@C 2711 PetscDSGetExactSolutionTimeDerivative - Get the pointwise time derivative of the exact solution function for a given test field 2712 2713 Not collective 2714 2715 Input Parameters: 2716 + prob - The `PetscDS` 2717 - f - The test field number 2718 2719 Output Parameters: 2720 + exactSol - time derivative of the exact solution for the test field 2721 - exactCtx - time derivative of the exact solution context 2722 2723 Calling sequence for the solution functions: 2724 .vb 2725 sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) 2726 .ve 2727 + dim - the spatial dimension 2728 . t - current time 2729 . x - coordinates of the current point 2730 . Nc - the number of field components 2731 . u - the solution field evaluated at the current point 2732 - ctx - a user context 2733 2734 Level: intermediate 2735 2736 .seealso: `PetscDS`, `PetscDSSetExactSolutionTimeDerivative()`, `PetscDSGetExactSolution()` 2737 @*/ 2738 PetscErrorCode PetscDSGetExactSolutionTimeDerivative(PetscDS prob, PetscInt f, PetscErrorCode (**sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void **ctx) 2739 { 2740 PetscFunctionBegin; 2741 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2742 PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf); 2743 if (sol) { 2744 PetscValidPointer(sol, 3); 2745 *sol = prob->exactSol_t[f]; 2746 } 2747 if (ctx) { 2748 PetscValidPointer(ctx, 4); 2749 *ctx = prob->exactCtx_t[f]; 2750 } 2751 PetscFunctionReturn(PETSC_SUCCESS); 2752 } 2753 2754 /*@C 2755 PetscDSSetExactSolutionTimeDerivative - Set the pointwise time derivative of the exact solution function for a given test field 2756 2757 Not collective 2758 2759 Input Parameters: 2760 + prob - The `PetscDS` 2761 . f - The test field number 2762 . sol - time derivative of the solution function for the test fields 2763 - ctx - time derivative of the solution context or NULL 2764 2765 Calling sequence for solution functions: 2766 .vb 2767 sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) 2768 .ve 2769 + dim - the spatial dimension 2770 . t - current time 2771 . x - coordinates of the current point 2772 . Nc - the number of field components 2773 . u - the solution field evaluated at the current point 2774 - ctx - a user context 2775 2776 Level: intermediate 2777 2778 .seealso: `PetscDS`, `PetscDSGetExactSolutionTimeDerivative()`, `PetscDSSetExactSolution()` 2779 @*/ 2780 PetscErrorCode PetscDSSetExactSolutionTimeDerivative(PetscDS prob, PetscInt f, PetscErrorCode (*sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void *ctx) 2781 { 2782 PetscFunctionBegin; 2783 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2784 PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f); 2785 PetscCall(PetscDSEnlarge_Static(prob, f + 1)); 2786 if (sol) { 2787 PetscValidFunction(sol, 3); 2788 prob->exactSol_t[f] = sol; 2789 } 2790 if (ctx) { 2791 PetscValidFunction(ctx, 4); 2792 prob->exactCtx_t[f] = ctx; 2793 } 2794 PetscFunctionReturn(PETSC_SUCCESS); 2795 } 2796 2797 /*@C 2798 PetscDSGetConstants - Returns the array of constants passed to point functions 2799 2800 Not collective 2801 2802 Input Parameter: 2803 . prob - The `PetscDS` object 2804 2805 Output Parameters: 2806 + numConstants - The number of constants 2807 - constants - The array of constants, NULL if there are none 2808 2809 Level: intermediate 2810 2811 .seealso: `PetscDS`, `PetscDSSetConstants()`, `PetscDSCreate()` 2812 @*/ 2813 PetscErrorCode PetscDSGetConstants(PetscDS prob, PetscInt *numConstants, const PetscScalar *constants[]) 2814 { 2815 PetscFunctionBegin; 2816 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2817 if (numConstants) { 2818 PetscValidIntPointer(numConstants, 2); 2819 *numConstants = prob->numConstants; 2820 } 2821 if (constants) { 2822 PetscValidPointer(constants, 3); 2823 *constants = prob->constants; 2824 } 2825 PetscFunctionReturn(PETSC_SUCCESS); 2826 } 2827 2828 /*@C 2829 PetscDSSetConstants - Set the array of constants passed to point functions 2830 2831 Not collective 2832 2833 Input Parameters: 2834 + prob - The `PetscDS` object 2835 . numConstants - The number of constants 2836 - constants - The array of constants, NULL if there are none 2837 2838 Level: intermediate 2839 2840 .seealso: `PetscDS`, `PetscDSGetConstants()`, `PetscDSCreate()` 2841 @*/ 2842 PetscErrorCode PetscDSSetConstants(PetscDS prob, PetscInt numConstants, PetscScalar constants[]) 2843 { 2844 PetscFunctionBegin; 2845 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2846 if (numConstants != prob->numConstants) { 2847 PetscCall(PetscFree(prob->constants)); 2848 prob->numConstants = numConstants; 2849 if (prob->numConstants) { 2850 PetscCall(PetscMalloc1(prob->numConstants, &prob->constants)); 2851 } else { 2852 prob->constants = NULL; 2853 } 2854 } 2855 if (prob->numConstants) { 2856 PetscValidScalarPointer(constants, 3); 2857 PetscCall(PetscArraycpy(prob->constants, constants, prob->numConstants)); 2858 } 2859 PetscFunctionReturn(PETSC_SUCCESS); 2860 } 2861 2862 /*@ 2863 PetscDSGetFieldIndex - Returns the index of the given field 2864 2865 Not collective 2866 2867 Input Parameters: 2868 + prob - The `PetscDS` object 2869 - disc - The discretization object 2870 2871 Output Parameter: 2872 . f - The field number 2873 2874 Level: beginner 2875 2876 .seealso: `PetscDS`, `PetscGetDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 2877 @*/ 2878 PetscErrorCode PetscDSGetFieldIndex(PetscDS prob, PetscObject disc, PetscInt *f) 2879 { 2880 PetscInt g; 2881 2882 PetscFunctionBegin; 2883 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2884 PetscValidIntPointer(f, 3); 2885 *f = -1; 2886 for (g = 0; g < prob->Nf; ++g) { 2887 if (disc == prob->disc[g]) break; 2888 } 2889 PetscCheck(g != prob->Nf, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Field not found in PetscDS."); 2890 *f = g; 2891 PetscFunctionReturn(PETSC_SUCCESS); 2892 } 2893 2894 /*@ 2895 PetscDSGetFieldSize - Returns the size of the given field in the full space basis 2896 2897 Not collective 2898 2899 Input Parameters: 2900 + prob - The `PetscDS` object 2901 - f - The field number 2902 2903 Output Parameter: 2904 . size - The size 2905 2906 Level: beginner 2907 2908 .seealso: `PetscDS`, `PetscDSGetFieldOffset()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 2909 @*/ 2910 PetscErrorCode PetscDSGetFieldSize(PetscDS prob, PetscInt f, PetscInt *size) 2911 { 2912 PetscFunctionBegin; 2913 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2914 PetscValidIntPointer(size, 3); 2915 PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf); 2916 PetscCall(PetscDSSetUp(prob)); 2917 *size = prob->Nb[f]; 2918 PetscFunctionReturn(PETSC_SUCCESS); 2919 } 2920 2921 /*@ 2922 PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis 2923 2924 Not collective 2925 2926 Input Parameters: 2927 + prob - The `PetscDS` object 2928 - f - The field number 2929 2930 Output Parameter: 2931 . off - The offset 2932 2933 Level: beginner 2934 2935 .seealso: `PetscDS`, `PetscDSGetFieldSize()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 2936 @*/ 2937 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off) 2938 { 2939 PetscInt size, g; 2940 2941 PetscFunctionBegin; 2942 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2943 PetscValidIntPointer(off, 3); 2944 PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf); 2945 *off = 0; 2946 for (g = 0; g < f; ++g) { 2947 PetscCall(PetscDSGetFieldSize(prob, g, &size)); 2948 *off += size; 2949 } 2950 PetscFunctionReturn(PETSC_SUCCESS); 2951 } 2952 2953 /*@ 2954 PetscDSGetFieldOffsetCohesive - Returns the offset of the given field in the full space basis on a cohesive cell 2955 2956 Not collective 2957 2958 Input Parameters: 2959 + prob - The `PetscDS` object 2960 - f - The field number 2961 2962 Output Parameter: 2963 . off - The offset 2964 2965 Level: beginner 2966 2967 .seealso: `PetscDS`, `PetscDSGetFieldSize()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 2968 @*/ 2969 PetscErrorCode PetscDSGetFieldOffsetCohesive(PetscDS ds, PetscInt f, PetscInt *off) 2970 { 2971 PetscInt size, g; 2972 2973 PetscFunctionBegin; 2974 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2975 PetscValidIntPointer(off, 3); 2976 PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf); 2977 *off = 0; 2978 for (g = 0; g < f; ++g) { 2979 PetscBool cohesive; 2980 2981 PetscCall(PetscDSGetCohesive(ds, g, &cohesive)); 2982 PetscCall(PetscDSGetFieldSize(ds, g, &size)); 2983 *off += cohesive ? size : size * 2; 2984 } 2985 PetscFunctionReturn(PETSC_SUCCESS); 2986 } 2987 2988 /*@ 2989 PetscDSGetDimensions - Returns the size of the approximation space for each field on an evaluation point 2990 2991 Not collective 2992 2993 Input Parameter: 2994 . prob - The `PetscDS` object 2995 2996 Output Parameter: 2997 . dimensions - The number of dimensions 2998 2999 Level: beginner 3000 3001 .seealso: `PetscDS`, `PetscDSGetComponentOffsets()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 3002 @*/ 3003 PetscErrorCode PetscDSGetDimensions(PetscDS prob, PetscInt *dimensions[]) 3004 { 3005 PetscFunctionBegin; 3006 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3007 PetscCall(PetscDSSetUp(prob)); 3008 PetscValidPointer(dimensions, 2); 3009 *dimensions = prob->Nb; 3010 PetscFunctionReturn(PETSC_SUCCESS); 3011 } 3012 3013 /*@ 3014 PetscDSGetComponents - Returns the number of components for each field on an evaluation point 3015 3016 Not collective 3017 3018 Input Parameter: 3019 . prob - The `PetscDS` object 3020 3021 Output Parameter: 3022 . components - The number of components 3023 3024 Level: beginner 3025 3026 .seealso: `PetscDS`, `PetscDSGetComponentOffsets()`, `PetscDSGetNumFields()`, `PetscDSCreate()` 3027 @*/ 3028 PetscErrorCode PetscDSGetComponents(PetscDS prob, PetscInt *components[]) 3029 { 3030 PetscFunctionBegin; 3031 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3032 PetscCall(PetscDSSetUp(prob)); 3033 PetscValidPointer(components, 2); 3034 *components = prob->Nc; 3035 PetscFunctionReturn(PETSC_SUCCESS); 3036 } 3037 3038 /*@ 3039 PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point 3040 3041 Not collective 3042 3043 Input Parameters: 3044 + prob - The `PetscDS` object 3045 - f - The field number 3046 3047 Output Parameter: 3048 . off - The offset 3049 3050 Level: beginner 3051 3052 .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()` 3053 @*/ 3054 PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off) 3055 { 3056 PetscFunctionBegin; 3057 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3058 PetscValidIntPointer(off, 3); 3059 PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf); 3060 PetscCall(PetscDSSetUp(prob)); 3061 *off = prob->off[f]; 3062 PetscFunctionReturn(PETSC_SUCCESS); 3063 } 3064 3065 /*@ 3066 PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point 3067 3068 Not collective 3069 3070 Input Parameter: 3071 . prob - The `PetscDS` object 3072 3073 Output Parameter: 3074 . offsets - The offsets 3075 3076 Level: beginner 3077 3078 .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()` 3079 @*/ 3080 PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[]) 3081 { 3082 PetscFunctionBegin; 3083 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3084 PetscValidPointer(offsets, 2); 3085 PetscCall(PetscDSSetUp(prob)); 3086 *offsets = prob->off; 3087 PetscFunctionReturn(PETSC_SUCCESS); 3088 } 3089 3090 /*@ 3091 PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point 3092 3093 Not collective 3094 3095 Input Parameter: 3096 . prob - The `PetscDS` object 3097 3098 Output Parameter: 3099 . offsets - The offsets 3100 3101 Level: beginner 3102 3103 .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()` 3104 @*/ 3105 PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[]) 3106 { 3107 PetscFunctionBegin; 3108 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3109 PetscValidPointer(offsets, 2); 3110 PetscCall(PetscDSSetUp(prob)); 3111 *offsets = prob->offDer; 3112 PetscFunctionReturn(PETSC_SUCCESS); 3113 } 3114 3115 /*@ 3116 PetscDSGetComponentOffsetsCohesive - Returns the offset of each field on an evaluation point 3117 3118 Not collective 3119 3120 Input Parameters: 3121 + ds - The `PetscDS` object 3122 - s - The cohesive side, 0 for negative, 1 for positive, 2 for cohesive 3123 3124 Output Parameter: 3125 . offsets - The offsets 3126 3127 Level: beginner 3128 3129 .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()` 3130 @*/ 3131 PetscErrorCode PetscDSGetComponentOffsetsCohesive(PetscDS ds, PetscInt s, PetscInt *offsets[]) 3132 { 3133 PetscFunctionBegin; 3134 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3135 PetscValidPointer(offsets, 3); 3136 PetscCheck(ds->isCohesive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cohesive offsets are only valid for a cohesive DS"); 3137 PetscCheck(!(s < 0) && !(s > 2), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cohesive side %" PetscInt_FMT " is not in [0, 2]", s); 3138 PetscCall(PetscDSSetUp(ds)); 3139 *offsets = ds->offCohesive[s]; 3140 PetscFunctionReturn(PETSC_SUCCESS); 3141 } 3142 3143 /*@ 3144 PetscDSGetComponentDerivativeOffsetsCohesive - Returns the offset of each field derivative on an evaluation point 3145 3146 Not collective 3147 3148 Input Parameters: 3149 + ds - The `PetscDS` object 3150 - s - The cohesive side, 0 for negative, 1 for positive, 2 for cohesive 3151 3152 Output Parameter: 3153 . offsets - The offsets 3154 3155 Level: beginner 3156 3157 .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()` 3158 @*/ 3159 PetscErrorCode PetscDSGetComponentDerivativeOffsetsCohesive(PetscDS ds, PetscInt s, PetscInt *offsets[]) 3160 { 3161 PetscFunctionBegin; 3162 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3163 PetscValidPointer(offsets, 3); 3164 PetscCheck(ds->isCohesive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cohesive offsets are only valid for a cohesive DS"); 3165 PetscCheck(!(s < 0) && !(s > 2), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cohesive side %" PetscInt_FMT " is not in [0, 2]", s); 3166 PetscCall(PetscDSSetUp(ds)); 3167 *offsets = ds->offDerCohesive[s]; 3168 PetscFunctionReturn(PETSC_SUCCESS); 3169 } 3170 3171 /*@C 3172 PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization 3173 3174 Not collective 3175 3176 Input Parameter: 3177 . prob - The `PetscDS` object 3178 3179 Output Parameter: 3180 . T - The basis function and derivatives tabulation at quadrature points for each field 3181 3182 Level: intermediate 3183 3184 .seealso: `PetscDS`, `PetscTabulation`, `PetscDSCreate()` 3185 @*/ 3186 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscTabulation *T[]) 3187 { 3188 PetscFunctionBegin; 3189 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3190 PetscValidPointer(T, 2); 3191 PetscCall(PetscDSSetUp(prob)); 3192 *T = prob->T; 3193 PetscFunctionReturn(PETSC_SUCCESS); 3194 } 3195 3196 /*@C 3197 PetscDSGetFaceTabulation - Return the basis tabulation at quadrature points on the faces 3198 3199 Not collective 3200 3201 Input Parameter: 3202 . prob - The `PetscDS` object 3203 3204 Output Parameter: 3205 . Tf - The basis function and derivative tabulation on each local face at quadrature points for each and field 3206 3207 Level: intermediate 3208 3209 .seealso: `PetscTabulation`, `PetscDS`, `PetscDSGetTabulation()`, `PetscDSCreate()` 3210 @*/ 3211 PetscErrorCode PetscDSGetFaceTabulation(PetscDS prob, PetscTabulation *Tf[]) 3212 { 3213 PetscFunctionBegin; 3214 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3215 PetscValidPointer(Tf, 2); 3216 PetscCall(PetscDSSetUp(prob)); 3217 *Tf = prob->Tf; 3218 PetscFunctionReturn(PETSC_SUCCESS); 3219 } 3220 3221 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x) 3222 { 3223 PetscFunctionBegin; 3224 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3225 PetscCall(PetscDSSetUp(prob)); 3226 if (u) { 3227 PetscValidPointer(u, 2); 3228 *u = prob->u; 3229 } 3230 if (u_t) { 3231 PetscValidPointer(u_t, 3); 3232 *u_t = prob->u_t; 3233 } 3234 if (u_x) { 3235 PetscValidPointer(u_x, 4); 3236 *u_x = prob->u_x; 3237 } 3238 PetscFunctionReturn(PETSC_SUCCESS); 3239 } 3240 3241 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3) 3242 { 3243 PetscFunctionBegin; 3244 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3245 PetscCall(PetscDSSetUp(prob)); 3246 if (f0) { 3247 PetscValidPointer(f0, 2); 3248 *f0 = prob->f0; 3249 } 3250 if (f1) { 3251 PetscValidPointer(f1, 3); 3252 *f1 = prob->f1; 3253 } 3254 if (g0) { 3255 PetscValidPointer(g0, 4); 3256 *g0 = prob->g0; 3257 } 3258 if (g1) { 3259 PetscValidPointer(g1, 5); 3260 *g1 = prob->g1; 3261 } 3262 if (g2) { 3263 PetscValidPointer(g2, 6); 3264 *g2 = prob->g2; 3265 } 3266 if (g3) { 3267 PetscValidPointer(g3, 7); 3268 *g3 = prob->g3; 3269 } 3270 PetscFunctionReturn(PETSC_SUCCESS); 3271 } 3272 3273 PetscErrorCode PetscDSGetWorkspace(PetscDS prob, PetscReal **x, PetscScalar **basisReal, PetscScalar **basisDerReal, PetscScalar **testReal, PetscScalar **testDerReal) 3274 { 3275 PetscFunctionBegin; 3276 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3277 PetscCall(PetscDSSetUp(prob)); 3278 if (x) { 3279 PetscValidPointer(x, 2); 3280 *x = prob->x; 3281 } 3282 if (basisReal) { 3283 PetscValidPointer(basisReal, 3); 3284 *basisReal = prob->basisReal; 3285 } 3286 if (basisDerReal) { 3287 PetscValidPointer(basisDerReal, 4); 3288 *basisDerReal = prob->basisDerReal; 3289 } 3290 if (testReal) { 3291 PetscValidPointer(testReal, 5); 3292 *testReal = prob->testReal; 3293 } 3294 if (testDerReal) { 3295 PetscValidPointer(testDerReal, 6); 3296 *testDerReal = prob->testDerReal; 3297 } 3298 PetscFunctionReturn(PETSC_SUCCESS); 3299 } 3300 3301 /*@C 3302 PetscDSAddBoundary - Add a boundary condition to the model. The pointwise functions are used to provide boundary values for essential boundary conditions. 3303 In FEM, they are acting upon by dual basis functionals to generate FEM coefficients which are fixed. Natural boundary conditions signal to PETSc that boundary 3304 integrals should be performed, using the kernels from `PetscDSSetBdResidual()`. 3305 3306 Collective on ds 3307 3308 Input Parameters: 3309 + ds - The PetscDS object 3310 . type - The type of condition, e.g. `DM_BC_ESSENTIAL`/`DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann) 3311 . name - The BC name 3312 . label - The label defining constrained points 3313 . Nv - The number of `DMLabel` values for constrained points 3314 . values - An array of label values for constrained points 3315 . field - The field to constrain 3316 . Nc - The number of constrained field components (0 will constrain all fields) 3317 . comps - An array of constrained component numbers 3318 . bcFunc - A pointwise function giving boundary values 3319 . bcFunc_t - A pointwise function giving the time derivative of the boundary values, or NULL 3320 - ctx - An optional user context for bcFunc 3321 3322 Output Parameters: 3323 - bd - The boundary number 3324 3325 Options Database Keys: 3326 + -bc_<boundary name> <num> - Overrides the boundary ids 3327 - -bc_<boundary name>_comp <num> - Overrides the boundary components 3328 3329 Level: developer 3330 3331 Note: 3332 Both bcFunc abd bcFunc_t will depend on the boundary condition type. If the type if `DM_BC_ESSENTIAL`, Then the calling sequence is: 3333 3334 $ bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[]) 3335 3336 If the type is `DM_BC_ESSENTIAL_FIELD` or other _FIELD value, then the calling sequence is: 3337 .vb 3338 bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, 3339 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 3340 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 3341 PetscReal time, const PetscReal x[], PetscScalar bcval[]) 3342 .ve 3343 + dim - the spatial dimension 3344 . Nf - the number of fields 3345 . uOff - the offset into u[] and u_t[] for each field 3346 . uOff_x - the offset into u_x[] for each field 3347 . u - each field evaluated at the current point 3348 . u_t - the time derivative of each field evaluated at the current point 3349 . u_x - the gradient of each field evaluated at the current point 3350 . aOff - the offset into a[] and a_t[] for each auxiliary field 3351 . aOff_x - the offset into a_x[] for each auxiliary field 3352 . a - each auxiliary field evaluated at the current point 3353 . a_t - the time derivative of each auxiliary field evaluated at the current point 3354 . a_x - the gradient of auxiliary each field evaluated at the current point 3355 . t - current time 3356 . x - coordinates of the current point 3357 . numConstants - number of constant parameters 3358 . constants - constant parameters 3359 - bcval - output values at the current point 3360 3361 .seealso: `PetscDS`, `PetscWeakForm`, `DMLabel`, `DMBoundaryConditionType`, `PetscDSAddBoundaryByName()`, `PetscDSGetBoundary()`, `PetscDSSetResidual()`, `PetscDSSetBdResidual()` 3362 @*/ 3363 PetscErrorCode PetscDSAddBoundary(PetscDS ds, DMBoundaryConditionType type, const char name[], DMLabel label, PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx, PetscInt *bd) 3364 { 3365 DSBoundary head = ds->boundary, b; 3366 PetscInt n = 0; 3367 const char *lname; 3368 3369 PetscFunctionBegin; 3370 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3371 PetscValidLogicalCollectiveEnum(ds, type, 2); 3372 PetscValidCharPointer(name, 3); 3373 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4); 3374 PetscValidLogicalCollectiveInt(ds, Nv, 5); 3375 PetscValidLogicalCollectiveInt(ds, field, 7); 3376 PetscValidLogicalCollectiveInt(ds, Nc, 8); 3377 PetscCheck(field >= 0 && field < ds->Nf, PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", field, ds->Nf); 3378 if (Nc > 0) { 3379 PetscInt *fcomps; 3380 PetscInt c; 3381 3382 PetscCall(PetscDSGetComponents(ds, &fcomps)); 3383 PetscCheck(Nc <= fcomps[field], PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_OUTOFRANGE, "Number of constrained components %" PetscInt_FMT " > %" PetscInt_FMT " components for field %" PetscInt_FMT, Nc, fcomps[field], field); 3384 for (c = 0; c < Nc; ++c) { 3385 PetscCheck(comps[c] >= 0 && comps[c] < fcomps[field], PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_OUTOFRANGE, "Constrained component[%" PetscInt_FMT "] %" PetscInt_FMT " not in [0, %" PetscInt_FMT ") components for field %" PetscInt_FMT, c, comps[c], fcomps[field], field); 3386 } 3387 } 3388 PetscCall(PetscNew(&b)); 3389 PetscCall(PetscStrallocpy(name, (char **)&b->name)); 3390 PetscCall(PetscWeakFormCreate(PETSC_COMM_SELF, &b->wf)); 3391 PetscCall(PetscWeakFormSetNumFields(b->wf, ds->Nf)); 3392 PetscCall(PetscMalloc1(Nv, &b->values)); 3393 if (Nv) PetscCall(PetscArraycpy(b->values, values, Nv)); 3394 PetscCall(PetscMalloc1(Nc, &b->comps)); 3395 if (Nc) PetscCall(PetscArraycpy(b->comps, comps, Nc)); 3396 PetscCall(PetscObjectGetName((PetscObject)label, &lname)); 3397 PetscCall(PetscStrallocpy(lname, (char **)&b->lname)); 3398 b->type = type; 3399 b->label = label; 3400 b->Nv = Nv; 3401 b->field = field; 3402 b->Nc = Nc; 3403 b->func = bcFunc; 3404 b->func_t = bcFunc_t; 3405 b->ctx = ctx; 3406 b->next = NULL; 3407 /* Append to linked list so that we can preserve the order */ 3408 if (!head) ds->boundary = b; 3409 while (head) { 3410 if (!head->next) { 3411 head->next = b; 3412 head = b; 3413 } 3414 head = head->next; 3415 ++n; 3416 } 3417 if (bd) { 3418 PetscValidIntPointer(bd, 13); 3419 *bd = n; 3420 } 3421 PetscFunctionReturn(PETSC_SUCCESS); 3422 } 3423 3424 /*@C 3425 PetscDSAddBoundaryByName - Add a boundary condition to the model. The pointwise functions are used to provide boundary values for essential boundary conditions. 3426 In FEM, they are acting upon by dual basis functionals to generate FEM coefficients which are fixed. Natural boundary conditions signal to PETSc that 3427 boundary integrals should be performed, using the kernels from `PetscDSSetBdResidual()`. 3428 3429 Collective on ds 3430 3431 Input Parameters: 3432 + ds - The `PetscDS` object 3433 . type - The type of condition, e.g. `DM_BC_ESSENTIAL`/`DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann) 3434 . name - The BC name 3435 . lname - The naem of the label defining constrained points 3436 . Nv - The number of `DMLabel` values for constrained points 3437 . values - An array of label values for constrained points 3438 . field - The field to constrain 3439 . Nc - The number of constrained field components (0 will constrain all fields) 3440 . comps - An array of constrained component numbers 3441 . bcFunc - A pointwise function giving boundary values 3442 . bcFunc_t - A pointwise function giving the time derivative of the boundary values, or NULL 3443 - ctx - An optional user context for bcFunc 3444 3445 Output Parameters: 3446 - bd - The boundary number 3447 3448 Options Database Keys: 3449 + -bc_<boundary name> <num> - Overrides the boundary ids 3450 - -bc_<boundary name>_comp <num> - Overrides the boundary components 3451 3452 Calling Sequence of bcFunc() and bcFunc_t(): 3453 If the type is `DM_BC_ESSENTIAL` 3454 .vb 3455 bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[]) 3456 .ve 3457 If the type is `DM_BC_ESSENTIAL_FIELD` or other _FIELD value, 3458 .vb 3459 bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, 3460 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 3461 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 3462 PetscReal time, const PetscReal x[], PetscScalar bcval[]) 3463 .ve 3464 + dim - the spatial dimension 3465 . Nf - the number of fields 3466 . uOff - the offset into u[] and u_t[] for each field 3467 . uOff_x - the offset into u_x[] for each field 3468 . u - each field evaluated at the current point 3469 . u_t - the time derivative of each field evaluated at the current point 3470 . u_x - the gradient of each field evaluated at the current point 3471 . aOff - the offset into a[] and a_t[] for each auxiliary field 3472 . aOff_x - the offset into a_x[] for each auxiliary field 3473 . a - each auxiliary field evaluated at the current point 3474 . a_t - the time derivative of each auxiliary field evaluated at the current point 3475 . a_x - the gradient of auxiliary each field evaluated at the current point 3476 . t - current time 3477 . x - coordinates of the current point 3478 . numConstants - number of constant parameters 3479 . constants - constant parameters 3480 - bcval - output values at the current point 3481 3482 Level: developer 3483 3484 Note: 3485 This function should only be used with `DMFOREST` currently, since labels cannot be defined before the underlying `DMPLEX` is built. 3486 3487 .seealso: `PetscDS`, `PetscWeakForm`, `DMLabel`, `DMBoundaryConditionType`, `PetscDSAddBoundary()`, `PetscDSGetBoundary()`, `PetscDSSetResidual()`, `PetscDSSetBdResidual()` 3488 @*/ 3489 PetscErrorCode PetscDSAddBoundaryByName(PetscDS ds, DMBoundaryConditionType type, const char name[], const char lname[], PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx, PetscInt *bd) 3490 { 3491 DSBoundary head = ds->boundary, b; 3492 PetscInt n = 0; 3493 3494 PetscFunctionBegin; 3495 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3496 PetscValidLogicalCollectiveEnum(ds, type, 2); 3497 PetscValidCharPointer(name, 3); 3498 PetscValidCharPointer(lname, 4); 3499 PetscValidLogicalCollectiveInt(ds, Nv, 5); 3500 PetscValidLogicalCollectiveInt(ds, field, 7); 3501 PetscValidLogicalCollectiveInt(ds, Nc, 8); 3502 PetscCall(PetscNew(&b)); 3503 PetscCall(PetscStrallocpy(name, (char **)&b->name)); 3504 PetscCall(PetscWeakFormCreate(PETSC_COMM_SELF, &b->wf)); 3505 PetscCall(PetscWeakFormSetNumFields(b->wf, ds->Nf)); 3506 PetscCall(PetscMalloc1(Nv, &b->values)); 3507 if (Nv) PetscCall(PetscArraycpy(b->values, values, Nv)); 3508 PetscCall(PetscMalloc1(Nc, &b->comps)); 3509 if (Nc) PetscCall(PetscArraycpy(b->comps, comps, Nc)); 3510 PetscCall(PetscStrallocpy(lname, (char **)&b->lname)); 3511 b->type = type; 3512 b->label = NULL; 3513 b->Nv = Nv; 3514 b->field = field; 3515 b->Nc = Nc; 3516 b->func = bcFunc; 3517 b->func_t = bcFunc_t; 3518 b->ctx = ctx; 3519 b->next = NULL; 3520 /* Append to linked list so that we can preserve the order */ 3521 if (!head) ds->boundary = b; 3522 while (head) { 3523 if (!head->next) { 3524 head->next = b; 3525 head = b; 3526 } 3527 head = head->next; 3528 ++n; 3529 } 3530 if (bd) { 3531 PetscValidIntPointer(bd, 13); 3532 *bd = n; 3533 } 3534 PetscFunctionReturn(PETSC_SUCCESS); 3535 } 3536 3537 /*@C 3538 PetscDSUpdateBoundary - Change a boundary condition for the model. The pointwise functions are used to provide boundary values for essential boundary conditions. 3539 In FEM, they are acting upon by dual basis functionals to generate FEM coefficients which are fixed. Natural boundary conditions signal to PETSc that boundary integrals 3540 should be performed, using the kernels from `PetscDSSetBdResidual()`. 3541 3542 Input Parameters: 3543 + ds - The `PetscDS` object 3544 . bd - The boundary condition number 3545 . type - The type of condition, e.g. `DM_BC_ESSENTIAL`/`DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann) 3546 . name - The BC name 3547 . label - The label defining constrained points 3548 . Nv - The number of `DMLabel` ids for constrained points 3549 . values - An array of ids for constrained points 3550 . field - The field to constrain 3551 . Nc - The number of constrained field components 3552 . comps - An array of constrained component numbers 3553 . bcFunc - A pointwise function giving boundary values 3554 . bcFunc_t - A pointwise function giving the time derivative of the boundary values, or NULL 3555 - ctx - An optional user context for bcFunc 3556 3557 Level: developer 3558 3559 Note: 3560 The boundary condition number is the order in which it was registered. The user can get the number of boundary conditions from `PetscDSGetNumBoundary()`. 3561 See `PetscDSAddBoundary()` for a description of the calling sequences for the callbacks. 3562 3563 .seealso: `PetscDS`, `PetscWeakForm`, `DMBoundaryConditionType`, `PetscDSAddBoundary()`, `PetscDSGetBoundary()`, `PetscDSGetNumBoundary()`, `DMLabel` 3564 @*/ 3565 PetscErrorCode PetscDSUpdateBoundary(PetscDS ds, PetscInt bd, DMBoundaryConditionType type, const char name[], DMLabel label, PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx) 3566 { 3567 DSBoundary b = ds->boundary; 3568 PetscInt n = 0; 3569 3570 PetscFunctionBegin; 3571 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3572 while (b) { 3573 if (n == bd) break; 3574 b = b->next; 3575 ++n; 3576 } 3577 PetscCheck(b, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", bd, n); 3578 if (name) { 3579 PetscCall(PetscFree(b->name)); 3580 PetscCall(PetscStrallocpy(name, (char **)&b->name)); 3581 } 3582 b->type = type; 3583 if (label) { 3584 const char *name; 3585 3586 b->label = label; 3587 PetscCall(PetscFree(b->lname)); 3588 PetscCall(PetscObjectGetName((PetscObject)label, &name)); 3589 PetscCall(PetscStrallocpy(name, (char **)&b->lname)); 3590 } 3591 if (Nv >= 0) { 3592 b->Nv = Nv; 3593 PetscCall(PetscFree(b->values)); 3594 PetscCall(PetscMalloc1(Nv, &b->values)); 3595 if (Nv) PetscCall(PetscArraycpy(b->values, values, Nv)); 3596 } 3597 if (field >= 0) b->field = field; 3598 if (Nc >= 0) { 3599 b->Nc = Nc; 3600 PetscCall(PetscFree(b->comps)); 3601 PetscCall(PetscMalloc1(Nc, &b->comps)); 3602 if (Nc) PetscCall(PetscArraycpy(b->comps, comps, Nc)); 3603 } 3604 if (bcFunc) b->func = bcFunc; 3605 if (bcFunc_t) b->func_t = bcFunc_t; 3606 if (ctx) b->ctx = ctx; 3607 PetscFunctionReturn(PETSC_SUCCESS); 3608 } 3609 3610 /*@ 3611 PetscDSGetNumBoundary - Get the number of registered BC 3612 3613 Input Parameters: 3614 . ds - The `PetscDS` object 3615 3616 Output Parameters: 3617 . numBd - The number of BC 3618 3619 Level: intermediate 3620 3621 .seealso: `PetscDS`, `PetscDSAddBoundary()`, `PetscDSGetBoundary()` 3622 @*/ 3623 PetscErrorCode PetscDSGetNumBoundary(PetscDS ds, PetscInt *numBd) 3624 { 3625 DSBoundary b = ds->boundary; 3626 3627 PetscFunctionBegin; 3628 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3629 PetscValidIntPointer(numBd, 2); 3630 *numBd = 0; 3631 while (b) { 3632 ++(*numBd); 3633 b = b->next; 3634 } 3635 PetscFunctionReturn(PETSC_SUCCESS); 3636 } 3637 3638 /*@C 3639 PetscDSGetBoundary - Gets a boundary condition to the model 3640 3641 Input Parameters: 3642 + ds - The `PetscDS` object 3643 - bd - The BC number 3644 3645 Output Parameters: 3646 + wf - The `PetscWeakForm` holding the pointwise functions 3647 . type - The type of condition, e.g. `DM_BC_ESSENTIAL`/`DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann) 3648 . name - The BC name 3649 . label - The label defining constrained points 3650 . Nv - The number of `DMLabel` ids for constrained points 3651 . values - An array of ids for constrained points 3652 . field - The field to constrain 3653 . Nc - The number of constrained field components 3654 . comps - An array of constrained component numbers 3655 . bcFunc - A pointwise function giving boundary values 3656 . bcFunc_t - A pointwise function giving the time derivative of the boundary values 3657 - ctx - An optional user context for bcFunc 3658 3659 Options Database Keys: 3660 + -bc_<boundary name> <num> - Overrides the boundary ids 3661 - -bc_<boundary name>_comp <num> - Overrides the boundary components 3662 3663 Level: developer 3664 3665 .seealso: `PetscDS`, `PetscWeakForm`, `DMBoundaryConditionType`, `PetscDSAddBoundary()`, `DMLabel` 3666 @*/ 3667 PetscErrorCode PetscDSGetBoundary(PetscDS ds, PetscInt bd, PetscWeakForm *wf, DMBoundaryConditionType *type, const char *name[], DMLabel *label, PetscInt *Nv, const PetscInt *values[], PetscInt *field, PetscInt *Nc, const PetscInt *comps[], void (**func)(void), void (**func_t)(void), void **ctx) 3668 { 3669 DSBoundary b = ds->boundary; 3670 PetscInt n = 0; 3671 3672 PetscFunctionBegin; 3673 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3674 while (b) { 3675 if (n == bd) break; 3676 b = b->next; 3677 ++n; 3678 } 3679 PetscCheck(b, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", bd, n); 3680 if (wf) { 3681 PetscValidPointer(wf, 3); 3682 *wf = b->wf; 3683 } 3684 if (type) { 3685 PetscValidPointer(type, 4); 3686 *type = b->type; 3687 } 3688 if (name) { 3689 PetscValidPointer(name, 5); 3690 *name = b->name; 3691 } 3692 if (label) { 3693 PetscValidPointer(label, 6); 3694 *label = b->label; 3695 } 3696 if (Nv) { 3697 PetscValidIntPointer(Nv, 7); 3698 *Nv = b->Nv; 3699 } 3700 if (values) { 3701 PetscValidPointer(values, 8); 3702 *values = b->values; 3703 } 3704 if (field) { 3705 PetscValidIntPointer(field, 9); 3706 *field = b->field; 3707 } 3708 if (Nc) { 3709 PetscValidIntPointer(Nc, 10); 3710 *Nc = b->Nc; 3711 } 3712 if (comps) { 3713 PetscValidPointer(comps, 11); 3714 *comps = b->comps; 3715 } 3716 if (func) { 3717 PetscValidPointer(func, 12); 3718 *func = b->func; 3719 } 3720 if (func_t) { 3721 PetscValidPointer(func_t, 13); 3722 *func_t = b->func_t; 3723 } 3724 if (ctx) { 3725 PetscValidPointer(ctx, 14); 3726 *ctx = b->ctx; 3727 } 3728 PetscFunctionReturn(PETSC_SUCCESS); 3729 } 3730 3731 static PetscErrorCode DSBoundaryDuplicate_Internal(DSBoundary b, DSBoundary *bNew) 3732 { 3733 PetscFunctionBegin; 3734 PetscCall(PetscNew(bNew)); 3735 PetscCall(PetscWeakFormCreate(PETSC_COMM_SELF, &(*bNew)->wf)); 3736 PetscCall(PetscWeakFormCopy(b->wf, (*bNew)->wf)); 3737 PetscCall(PetscStrallocpy(b->name, (char **)&((*bNew)->name))); 3738 PetscCall(PetscStrallocpy(b->lname, (char **)&((*bNew)->lname))); 3739 (*bNew)->type = b->type; 3740 (*bNew)->label = b->label; 3741 (*bNew)->Nv = b->Nv; 3742 PetscCall(PetscMalloc1(b->Nv, &(*bNew)->values)); 3743 PetscCall(PetscArraycpy((*bNew)->values, b->values, b->Nv)); 3744 (*bNew)->field = b->field; 3745 (*bNew)->Nc = b->Nc; 3746 PetscCall(PetscMalloc1(b->Nc, &(*bNew)->comps)); 3747 PetscCall(PetscArraycpy((*bNew)->comps, b->comps, b->Nc)); 3748 (*bNew)->func = b->func; 3749 (*bNew)->func_t = b->func_t; 3750 (*bNew)->ctx = b->ctx; 3751 PetscFunctionReturn(PETSC_SUCCESS); 3752 } 3753 3754 /*@ 3755 PetscDSCopyBoundary - Copy all boundary condition objects to the new problem 3756 3757 Not collective 3758 3759 Input Parameters: 3760 + ds - The source `PetscDS` object 3761 . numFields - The number of selected fields, or `PETSC_DEFAULT` for all fields 3762 - fields - The selected fields, or NULL for all fields 3763 3764 Output Parameter: 3765 . newds - The target `PetscDS`, now with a copy of the boundary conditions 3766 3767 Level: intermediate 3768 3769 .seealso: `PetscDS`, `DMBoundary`, `PetscDSCopyEquations()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()` 3770 @*/ 3771 PetscErrorCode PetscDSCopyBoundary(PetscDS ds, PetscInt numFields, const PetscInt fields[], PetscDS newds) 3772 { 3773 DSBoundary b, *lastnext; 3774 3775 PetscFunctionBegin; 3776 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3777 PetscValidHeaderSpecific(newds, PETSCDS_CLASSID, 4); 3778 if (ds == newds) PetscFunctionReturn(PETSC_SUCCESS); 3779 PetscCall(PetscDSDestroyBoundary(newds)); 3780 lastnext = &(newds->boundary); 3781 for (b = ds->boundary; b; b = b->next) { 3782 DSBoundary bNew; 3783 PetscInt fieldNew = -1; 3784 3785 if (numFields > 0 && fields) { 3786 PetscInt f; 3787 3788 for (f = 0; f < numFields; ++f) 3789 if (b->field == fields[f]) break; 3790 if (f == numFields) continue; 3791 fieldNew = f; 3792 } 3793 PetscCall(DSBoundaryDuplicate_Internal(b, &bNew)); 3794 bNew->field = fieldNew < 0 ? b->field : fieldNew; 3795 *lastnext = bNew; 3796 lastnext = &(bNew->next); 3797 } 3798 PetscFunctionReturn(PETSC_SUCCESS); 3799 } 3800 3801 /*@ 3802 PetscDSDestroyBoundary - Remove all `DMBoundary` objects from the `PetscDS` 3803 3804 Not collective 3805 3806 Input Parameter: 3807 . ds - The `PetscDS` object 3808 3809 Level: intermediate 3810 3811 .seealso: `PetscDS`, `DMBoundary`, `PetscDSCopyBoundary()`, `PetscDSCopyEquations()` 3812 @*/ 3813 PetscErrorCode PetscDSDestroyBoundary(PetscDS ds) 3814 { 3815 DSBoundary next = ds->boundary; 3816 3817 PetscFunctionBegin; 3818 while (next) { 3819 DSBoundary b = next; 3820 3821 next = b->next; 3822 PetscCall(PetscWeakFormDestroy(&b->wf)); 3823 PetscCall(PetscFree(b->name)); 3824 PetscCall(PetscFree(b->lname)); 3825 PetscCall(PetscFree(b->values)); 3826 PetscCall(PetscFree(b->comps)); 3827 PetscCall(PetscFree(b)); 3828 } 3829 PetscFunctionReturn(PETSC_SUCCESS); 3830 } 3831 3832 /*@ 3833 PetscDSSelectDiscretizations - Copy discretizations to the new problem with different field layout 3834 3835 Not collective 3836 3837 Input Parameters: 3838 + prob - The `PetscDS` object 3839 . numFields - Number of new fields 3840 - fields - Old field number for each new field 3841 3842 Output Parameter: 3843 . newprob - The `PetscDS` copy 3844 3845 Level: intermediate 3846 3847 .seealso: `PetscDS`, `PetscDSSelectEquations()`, `PetscDSCopyBoundary()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()` 3848 @*/ 3849 PetscErrorCode PetscDSSelectDiscretizations(PetscDS prob, PetscInt numFields, const PetscInt fields[], PetscDS newprob) 3850 { 3851 PetscInt Nf, Nfn, fn; 3852 3853 PetscFunctionBegin; 3854 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3855 if (fields) PetscValidIntPointer(fields, 3); 3856 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 4); 3857 PetscCall(PetscDSGetNumFields(prob, &Nf)); 3858 PetscCall(PetscDSGetNumFields(newprob, &Nfn)); 3859 numFields = numFields < 0 ? Nf : numFields; 3860 for (fn = 0; fn < numFields; ++fn) { 3861 const PetscInt f = fields ? fields[fn] : fn; 3862 PetscObject disc; 3863 3864 if (f >= Nf) continue; 3865 PetscCall(PetscDSGetDiscretization(prob, f, &disc)); 3866 PetscCall(PetscDSSetDiscretization(newprob, fn, disc)); 3867 } 3868 PetscFunctionReturn(PETSC_SUCCESS); 3869 } 3870 3871 /*@ 3872 PetscDSSelectEquations - Copy pointwise function pointers to the new problem with different field layout 3873 3874 Not collective 3875 3876 Input Parameters: 3877 + prob - The `PetscDS` object 3878 . numFields - Number of new fields 3879 - fields - Old field number for each new field 3880 3881 Output Parameter: 3882 . newprob - The `PetscDS` copy 3883 3884 Level: intermediate 3885 3886 .seealso: `PetscDS`, `PetscDSSelectDiscretizations()`, `PetscDSCopyBoundary()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()` 3887 @*/ 3888 PetscErrorCode PetscDSSelectEquations(PetscDS prob, PetscInt numFields, const PetscInt fields[], PetscDS newprob) 3889 { 3890 PetscInt Nf, Nfn, fn, gn; 3891 3892 PetscFunctionBegin; 3893 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3894 if (fields) PetscValidIntPointer(fields, 3); 3895 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 4); 3896 PetscCall(PetscDSGetNumFields(prob, &Nf)); 3897 PetscCall(PetscDSGetNumFields(newprob, &Nfn)); 3898 PetscCheck(numFields <= Nfn, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_SIZ, "Number of fields %" PetscInt_FMT " to transfer must not be greater then the total number of fields %" PetscInt_FMT, numFields, Nfn); 3899 for (fn = 0; fn < numFields; ++fn) { 3900 const PetscInt f = fields ? fields[fn] : fn; 3901 PetscPointFunc obj; 3902 PetscPointFunc f0, f1; 3903 PetscBdPointFunc f0Bd, f1Bd; 3904 PetscRiemannFunc r; 3905 3906 if (f >= Nf) continue; 3907 PetscCall(PetscDSGetObjective(prob, f, &obj)); 3908 PetscCall(PetscDSGetResidual(prob, f, &f0, &f1)); 3909 PetscCall(PetscDSGetBdResidual(prob, f, &f0Bd, &f1Bd)); 3910 PetscCall(PetscDSGetRiemannSolver(prob, f, &r)); 3911 PetscCall(PetscDSSetObjective(newprob, fn, obj)); 3912 PetscCall(PetscDSSetResidual(newprob, fn, f0, f1)); 3913 PetscCall(PetscDSSetBdResidual(newprob, fn, f0Bd, f1Bd)); 3914 PetscCall(PetscDSSetRiemannSolver(newprob, fn, r)); 3915 for (gn = 0; gn < numFields; ++gn) { 3916 const PetscInt g = fields ? fields[gn] : gn; 3917 PetscPointJac g0, g1, g2, g3; 3918 PetscPointJac g0p, g1p, g2p, g3p; 3919 PetscBdPointJac g0Bd, g1Bd, g2Bd, g3Bd; 3920 3921 if (g >= Nf) continue; 3922 PetscCall(PetscDSGetJacobian(prob, f, g, &g0, &g1, &g2, &g3)); 3923 PetscCall(PetscDSGetJacobianPreconditioner(prob, f, g, &g0p, &g1p, &g2p, &g3p)); 3924 PetscCall(PetscDSGetBdJacobian(prob, f, g, &g0Bd, &g1Bd, &g2Bd, &g3Bd)); 3925 PetscCall(PetscDSSetJacobian(newprob, fn, gn, g0, g1, g2, g3)); 3926 PetscCall(PetscDSSetJacobianPreconditioner(newprob, fn, gn, g0p, g1p, g2p, g3p)); 3927 PetscCall(PetscDSSetBdJacobian(newprob, fn, gn, g0Bd, g1Bd, g2Bd, g3Bd)); 3928 } 3929 } 3930 PetscFunctionReturn(PETSC_SUCCESS); 3931 } 3932 3933 /*@ 3934 PetscDSCopyEquations - Copy all pointwise function pointers to another `PetscDS` 3935 3936 Not collective 3937 3938 Input Parameter: 3939 . prob - The `PetscDS` object 3940 3941 Output Parameter: 3942 . newprob - The `PetscDS` copy 3943 3944 Level: intermediate 3945 3946 .seealso: `PetscDS`, `PetscDSCopyBoundary()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()` 3947 @*/ 3948 PetscErrorCode PetscDSCopyEquations(PetscDS prob, PetscDS newprob) 3949 { 3950 PetscWeakForm wf, newwf; 3951 PetscInt Nf, Ng; 3952 3953 PetscFunctionBegin; 3954 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3955 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2); 3956 PetscCall(PetscDSGetNumFields(prob, &Nf)); 3957 PetscCall(PetscDSGetNumFields(newprob, &Ng)); 3958 PetscCheck(Nf == Ng, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_SIZ, "Number of fields must match %" PetscInt_FMT " != %" PetscInt_FMT, Nf, Ng); 3959 PetscCall(PetscDSGetWeakForm(prob, &wf)); 3960 PetscCall(PetscDSGetWeakForm(newprob, &newwf)); 3961 PetscCall(PetscWeakFormCopy(wf, newwf)); 3962 PetscFunctionReturn(PETSC_SUCCESS); 3963 } 3964 3965 /*@ 3966 PetscDSCopyConstants - Copy all constants to another `PetscDS` 3967 3968 Not collective 3969 3970 Input Parameter: 3971 . prob - The `PetscDS` object 3972 3973 Output Parameter: 3974 . newprob - The `PetscDS` copy 3975 3976 Level: intermediate 3977 3978 .seealso: `PetscDS`, `PetscDSCopyBoundary()`, `PetscDSCopyEquations()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()` 3979 @*/ 3980 PetscErrorCode PetscDSCopyConstants(PetscDS prob, PetscDS newprob) 3981 { 3982 PetscInt Nc; 3983 const PetscScalar *constants; 3984 3985 PetscFunctionBegin; 3986 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3987 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2); 3988 PetscCall(PetscDSGetConstants(prob, &Nc, &constants)); 3989 PetscCall(PetscDSSetConstants(newprob, Nc, (PetscScalar *)constants)); 3990 PetscFunctionReturn(PETSC_SUCCESS); 3991 } 3992 3993 /*@ 3994 PetscDSCopyExactSolutions - Copy all exact solutions to another `PetscDS` 3995 3996 Not collective 3997 3998 Input Parameter: 3999 . ds - The `PetscDS` object 4000 4001 Output Parameter: 4002 . newds - The `PetscDS` copy 4003 4004 Level: intermediate 4005 4006 .seealso: `PetscDS`, `PetscDSCopyBoundary()`, `PetscDSCopyEquations()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()` 4007 @*/ 4008 PetscErrorCode PetscDSCopyExactSolutions(PetscDS ds, PetscDS newds) 4009 { 4010 PetscSimplePointFunc sol; 4011 void *ctx; 4012 PetscInt Nf, f; 4013 4014 PetscFunctionBegin; 4015 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 4016 PetscValidHeaderSpecific(newds, PETSCDS_CLASSID, 2); 4017 PetscCall(PetscDSGetNumFields(ds, &Nf)); 4018 for (f = 0; f < Nf; ++f) { 4019 PetscCall(PetscDSGetExactSolution(ds, f, &sol, &ctx)); 4020 PetscCall(PetscDSSetExactSolution(newds, f, sol, ctx)); 4021 PetscCall(PetscDSGetExactSolutionTimeDerivative(ds, f, &sol, &ctx)); 4022 PetscCall(PetscDSSetExactSolutionTimeDerivative(newds, f, sol, ctx)); 4023 } 4024 PetscFunctionReturn(PETSC_SUCCESS); 4025 } 4026 4027 PetscErrorCode PetscDSGetHeightSubspace(PetscDS prob, PetscInt height, PetscDS *subprob) 4028 { 4029 PetscInt dim, Nf, f; 4030 4031 PetscFunctionBegin; 4032 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 4033 PetscValidPointer(subprob, 3); 4034 if (height == 0) { 4035 *subprob = prob; 4036 PetscFunctionReturn(PETSC_SUCCESS); 4037 } 4038 PetscCall(PetscDSGetNumFields(prob, &Nf)); 4039 PetscCall(PetscDSGetSpatialDimension(prob, &dim)); 4040 PetscCheck(height <= dim, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_OUTOFRANGE, "DS can only handle height in [0, %" PetscInt_FMT "], not %" PetscInt_FMT, dim, height); 4041 if (!prob->subprobs) PetscCall(PetscCalloc1(dim, &prob->subprobs)); 4042 if (!prob->subprobs[height - 1]) { 4043 PetscInt cdim; 4044 4045 PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)prob), &prob->subprobs[height - 1])); 4046 PetscCall(PetscDSGetCoordinateDimension(prob, &cdim)); 4047 PetscCall(PetscDSSetCoordinateDimension(prob->subprobs[height - 1], cdim)); 4048 for (f = 0; f < Nf; ++f) { 4049 PetscFE subfe; 4050 PetscObject obj; 4051 PetscClassId id; 4052 4053 PetscCall(PetscDSGetDiscretization(prob, f, &obj)); 4054 PetscCall(PetscObjectGetClassId(obj, &id)); 4055 if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetHeightSubspace((PetscFE)obj, height, &subfe)); 4056 else SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unsupported discretization type for field %" PetscInt_FMT, f); 4057 PetscCall(PetscDSSetDiscretization(prob->subprobs[height - 1], f, (PetscObject)subfe)); 4058 } 4059 } 4060 *subprob = prob->subprobs[height - 1]; 4061 PetscFunctionReturn(PETSC_SUCCESS); 4062 } 4063 4064 PetscErrorCode PetscDSGetDiscType_Internal(PetscDS ds, PetscInt f, PetscDiscType *disctype) 4065 { 4066 PetscObject obj; 4067 PetscClassId id; 4068 PetscInt Nf; 4069 4070 PetscFunctionBegin; 4071 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 4072 PetscValidPointer(disctype, 3); 4073 *disctype = PETSC_DISC_NONE; 4074 PetscCall(PetscDSGetNumFields(ds, &Nf)); 4075 PetscCheck(f < Nf, PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_SIZ, "Field %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, Nf); 4076 PetscCall(PetscDSGetDiscretization(ds, f, &obj)); 4077 if (obj) { 4078 PetscCall(PetscObjectGetClassId(obj, &id)); 4079 if (id == PETSCFE_CLASSID) *disctype = PETSC_DISC_FE; 4080 else *disctype = PETSC_DISC_FV; 4081 } 4082 PetscFunctionReturn(PETSC_SUCCESS); 4083 } 4084 4085 static PetscErrorCode PetscDSDestroy_Basic(PetscDS ds) 4086 { 4087 PetscFunctionBegin; 4088 PetscCall(PetscFree(ds->data)); 4089 PetscFunctionReturn(PETSC_SUCCESS); 4090 } 4091 4092 static PetscErrorCode PetscDSInitialize_Basic(PetscDS ds) 4093 { 4094 PetscFunctionBegin; 4095 ds->ops->setfromoptions = NULL; 4096 ds->ops->setup = NULL; 4097 ds->ops->view = NULL; 4098 ds->ops->destroy = PetscDSDestroy_Basic; 4099 PetscFunctionReturn(PETSC_SUCCESS); 4100 } 4101 4102 /*MC 4103 PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions 4104 4105 Level: intermediate 4106 4107 .seealso: `PetscDSType`, `PetscDSCreate()`, `PetscDSSetType()` 4108 M*/ 4109 4110 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS ds) 4111 { 4112 PetscDS_Basic *b; 4113 4114 PetscFunctionBegin; 4115 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 4116 PetscCall(PetscNew(&b)); 4117 ds->data = b; 4118 4119 PetscCall(PetscDSInitialize_Basic(ds)); 4120 PetscFunctionReturn(PETSC_SUCCESS); 4121 } 4122