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 Notes: 31 PetscDSRegister() may be called multiple times to add several user-defined PetscDSs 32 33 Sample usage: 34 .vb 35 PetscDSRegister("my_ds", MyPetscDSCreate); 36 .ve 37 38 Then, your PetscDS type can be chosen with the procedural interface via 39 .vb 40 PetscDSCreate(MPI_Comm, PetscDS *); 41 PetscDSSetType(PetscDS, "my_ds"); 42 .ve 43 or at runtime via the option 44 .vb 45 -petscds_type my_ds 46 .ve 47 48 Level: advanced 49 50 Not available from Fortran 51 52 .keywords: PetscDS, register 53 .seealso: PetscDSRegisterAll(), PetscDSRegisterDestroy() 54 55 @*/ 56 PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS)) 57 { 58 PetscErrorCode ierr; 59 60 PetscFunctionBegin; 61 ierr = PetscFunctionListAdd(&PetscDSList, sname, function);CHKERRQ(ierr); 62 PetscFunctionReturn(0); 63 } 64 65 /*@C 66 PetscDSSetType - Builds a particular PetscDS 67 68 Collective on PetscDS 69 70 Input Parameters: 71 + prob - The PetscDS object 72 - name - The kind of system 73 74 Options Database Key: 75 . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types 76 77 Level: intermediate 78 79 Not available from Fortran 80 81 .keywords: PetscDS, set, type 82 .seealso: PetscDSGetType(), PetscDSCreate() 83 @*/ 84 PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name) 85 { 86 PetscErrorCode (*r)(PetscDS); 87 PetscBool match; 88 PetscErrorCode ierr; 89 90 PetscFunctionBegin; 91 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 92 ierr = PetscObjectTypeCompare((PetscObject) prob, name, &match);CHKERRQ(ierr); 93 if (match) PetscFunctionReturn(0); 94 95 ierr = PetscDSRegisterAll();CHKERRQ(ierr); 96 ierr = PetscFunctionListFind(PetscDSList, name, &r);CHKERRQ(ierr); 97 if (!r) SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name); 98 99 if (prob->ops->destroy) { 100 ierr = (*prob->ops->destroy)(prob);CHKERRQ(ierr); 101 prob->ops->destroy = NULL; 102 } 103 ierr = (*r)(prob);CHKERRQ(ierr); 104 ierr = PetscObjectChangeTypeName((PetscObject) prob, name);CHKERRQ(ierr); 105 PetscFunctionReturn(0); 106 } 107 108 /*@C 109 PetscDSGetType - Gets the PetscDS type name (as a string) from the object. 110 111 Not Collective 112 113 Input Parameter: 114 . prob - The PetscDS 115 116 Output Parameter: 117 . name - The PetscDS type name 118 119 Level: intermediate 120 121 Not available from Fortran 122 123 .keywords: PetscDS, get, type, name 124 .seealso: PetscDSSetType(), PetscDSCreate() 125 @*/ 126 PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name) 127 { 128 PetscErrorCode ierr; 129 130 PetscFunctionBegin; 131 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 132 PetscValidPointer(name, 2); 133 ierr = PetscDSRegisterAll();CHKERRQ(ierr); 134 *name = ((PetscObject) prob)->type_name; 135 PetscFunctionReturn(0); 136 } 137 138 static PetscErrorCode PetscDSView_Ascii(PetscDS prob, PetscViewer viewer) 139 { 140 PetscViewerFormat format; 141 const PetscScalar *constants; 142 PetscInt numConstants, f; 143 PetscErrorCode ierr; 144 145 PetscFunctionBegin; 146 ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 147 ierr = PetscViewerASCIIPrintf(viewer, "Discrete System with %d fields\n", prob->Nf);CHKERRQ(ierr); 148 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 149 ierr = PetscViewerASCIIPrintf(viewer, " cell total dim %D total comp %D\n", prob->totDim, prob->totComp);CHKERRQ(ierr); 150 if (prob->isHybrid) {ierr = PetscViewerASCIIPrintf(viewer, " hybrid cell\n");CHKERRQ(ierr);} 151 for (f = 0; f < prob->Nf; ++f) { 152 DSBoundary b; 153 PetscObject obj; 154 PetscClassId id; 155 PetscQuadrature q; 156 const char *name; 157 PetscInt Nc, Nq, Nqc; 158 159 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 160 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 161 ierr = PetscObjectGetName(obj, &name);CHKERRQ(ierr); 162 ierr = PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>");CHKERRQ(ierr); 163 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);CHKERRQ(ierr); 164 if (id == PETSCFE_CLASSID) { 165 ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr); 166 ierr = PetscFEGetQuadrature((PetscFE) obj, &q);CHKERRQ(ierr); 167 ierr = PetscViewerASCIIPrintf(viewer, " FEM");CHKERRQ(ierr); 168 } else if (id == PETSCFV_CLASSID) { 169 ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr); 170 ierr = PetscFVGetQuadrature((PetscFV) obj, &q);CHKERRQ(ierr); 171 ierr = PetscViewerASCIIPrintf(viewer, " FVM");CHKERRQ(ierr); 172 } 173 else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f); 174 if (Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "%D components", Nc);CHKERRQ(ierr);} 175 else {ierr = PetscViewerASCIIPrintf(viewer, "%D component ", Nc);CHKERRQ(ierr);} 176 if (prob->implicit[f]) {ierr = PetscViewerASCIIPrintf(viewer, " (implicit)");CHKERRQ(ierr);} 177 else {ierr = PetscViewerASCIIPrintf(viewer, " (explicit)");CHKERRQ(ierr);} 178 if (q) { 179 ierr = PetscQuadratureGetData(q, NULL, &Nqc, &Nq, NULL, NULL);CHKERRQ(ierr); 180 ierr = PetscViewerASCIIPrintf(viewer, " (Nq %D Nqc %D)", Nq, Nqc);CHKERRQ(ierr); 181 } 182 ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 183 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);CHKERRQ(ierr); 184 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 185 if (id == PETSCFE_CLASSID) {ierr = PetscFEView((PetscFE) obj, viewer);CHKERRQ(ierr);} 186 else if (id == PETSCFV_CLASSID) {ierr = PetscFVView((PetscFV) obj, viewer);CHKERRQ(ierr);} 187 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 188 189 for (b = prob->boundary; b; b = b->next) { 190 PetscInt c, i; 191 192 if (b->field != f) continue; 193 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 194 ierr = PetscViewerASCIIPrintf(viewer, "Boundary %s (%s) %s\n", b->name, b->labelname, DMBoundaryConditionTypes[b->type]);CHKERRQ(ierr); 195 if (!b->numcomps) { 196 ierr = PetscViewerASCIIPrintf(viewer, " all components\n");CHKERRQ(ierr); 197 } else { 198 ierr = PetscViewerASCIIPrintf(viewer, " components: ");CHKERRQ(ierr); 199 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);CHKERRQ(ierr); 200 for (c = 0; c < b->numcomps; ++c) { 201 if (c > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 202 ierr = PetscViewerASCIIPrintf(viewer, "%D", b->comps[c]);CHKERRQ(ierr); 203 } 204 ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 205 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);CHKERRQ(ierr); 206 } 207 ierr = PetscViewerASCIIPrintf(viewer, " ids: ");CHKERRQ(ierr); 208 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);CHKERRQ(ierr); 209 for (i = 0; i < b->numids; ++i) { 210 if (i > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 211 ierr = PetscViewerASCIIPrintf(viewer, "%D", b->ids[i]);CHKERRQ(ierr); 212 } 213 ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 214 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);CHKERRQ(ierr); 215 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 216 } 217 } 218 ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 219 if (numConstants) { 220 ierr = PetscViewerASCIIPrintf(viewer, "%D constants\n", numConstants);CHKERRQ(ierr); 221 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 222 for (f = 0; f < numConstants; ++f) {ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double) PetscRealPart(constants[f]));CHKERRQ(ierr);} 223 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 224 } 225 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 226 PetscFunctionReturn(0); 227 } 228 229 /*@C 230 PetscDSView - Views a PetscDS 231 232 Collective on PetscDS 233 234 Input Parameter: 235 + prob - the PetscDS object to view 236 - v - the viewer 237 238 Level: developer 239 240 .seealso PetscDSDestroy() 241 @*/ 242 PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v) 243 { 244 PetscBool iascii; 245 PetscErrorCode ierr; 246 247 PetscFunctionBegin; 248 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 249 if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);} 250 else {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);} 251 ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 252 if (iascii) {ierr = PetscDSView_Ascii(prob, v);CHKERRQ(ierr);} 253 if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);} 254 PetscFunctionReturn(0); 255 } 256 257 /*@ 258 PetscDSSetFromOptions - sets parameters in a PetscDS from the options database 259 260 Collective on PetscDS 261 262 Input Parameter: 263 . prob - the PetscDS object to set options for 264 265 Options Database: 266 + -petscds_type <type> : Set the DS type 267 . -petscds_view <view opt> : View the DS 268 . -petscds_jac_pre : Turn formation of a separate Jacobian preconditioner on and off 269 . -bc_<name> <ids> : Specify a list of label ids for a boundary condition 270 - -bc_<name>_comp <comps> : Specify a list of field components to constrain for a boundary condition 271 272 Level: developer 273 274 .seealso PetscDSView() 275 @*/ 276 PetscErrorCode PetscDSSetFromOptions(PetscDS prob) 277 { 278 DSBoundary b; 279 const char *defaultType; 280 char name[256]; 281 PetscBool flg; 282 PetscErrorCode ierr; 283 284 PetscFunctionBegin; 285 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 286 if (!((PetscObject) prob)->type_name) { 287 defaultType = PETSCDSBASIC; 288 } else { 289 defaultType = ((PetscObject) prob)->type_name; 290 } 291 ierr = PetscDSRegisterAll();CHKERRQ(ierr); 292 293 ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr); 294 for (b = prob->boundary; b; b = b->next) { 295 char optname[1024]; 296 PetscInt ids[1024], len = 1024; 297 PetscBool flg; 298 299 ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr); 300 ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr); 301 ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr); 302 if (flg) { 303 b->numids = len; 304 ierr = PetscFree(b->ids);CHKERRQ(ierr); 305 ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr); 306 ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr); 307 } 308 len = 1024; 309 ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s_comp", b->name);CHKERRQ(ierr); 310 ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr); 311 ierr = PetscOptionsIntArray(optname, "List of boundary field components", "", ids, &len, &flg);CHKERRQ(ierr); 312 if (flg) { 313 b->numcomps = len; 314 ierr = PetscFree(b->comps);CHKERRQ(ierr); 315 ierr = PetscMalloc1(len, &b->comps);CHKERRQ(ierr); 316 ierr = PetscMemcpy(b->comps, ids, len*sizeof(PetscInt));CHKERRQ(ierr); 317 } 318 } 319 ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr); 320 if (flg) { 321 ierr = PetscDSSetType(prob, name);CHKERRQ(ierr); 322 } else if (!((PetscObject) prob)->type_name) { 323 ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr); 324 } 325 ierr = PetscOptionsBool("-petscds_jac_pre", "Discrete System", "PetscDSUseJacobianPreconditioner", prob->useJacPre, &prob->useJacPre, &flg);CHKERRQ(ierr); 326 if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);} 327 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 328 ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) prob);CHKERRQ(ierr); 329 ierr = PetscOptionsEnd();CHKERRQ(ierr); 330 if (prob->Nf) {ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr);} 331 PetscFunctionReturn(0); 332 } 333 334 /*@C 335 PetscDSSetUp - Construct data structures for the PetscDS 336 337 Collective on PetscDS 338 339 Input Parameter: 340 . prob - the PetscDS object to setup 341 342 Level: developer 343 344 .seealso PetscDSView(), PetscDSDestroy() 345 @*/ 346 PetscErrorCode PetscDSSetUp(PetscDS prob) 347 { 348 const PetscInt Nf = prob->Nf; 349 PetscInt dim, dimEmbed, NbMax = 0, NcMax = 0, NqMax = 0, NsMax = 1, f; 350 PetscErrorCode ierr; 351 352 PetscFunctionBegin; 353 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 354 if (prob->setup) PetscFunctionReturn(0); 355 /* Calculate sizes */ 356 ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr); 357 ierr = PetscDSGetCoordinateDimension(prob, &dimEmbed);CHKERRQ(ierr); 358 prob->totDim = prob->totComp = 0; 359 ierr = PetscMalloc2(Nf,&prob->Nc,Nf,&prob->Nb);CHKERRQ(ierr); 360 ierr = PetscCalloc2(Nf+1,&prob->off,Nf+1,&prob->offDer);CHKERRQ(ierr); 361 ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisFace,Nf,&prob->basisDerFace);CHKERRQ(ierr); 362 for (f = 0; f < Nf; ++f) { 363 PetscObject obj; 364 PetscClassId id; 365 PetscQuadrature q; 366 PetscInt Nq = 0, Nb, Nc; 367 368 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 369 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 370 if (id == PETSCFE_CLASSID) { 371 PetscFE fe = (PetscFE) obj; 372 373 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 374 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 375 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 376 ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr); 377 ierr = PetscFEGetFaceTabulation(fe, &prob->basisFace[f], &prob->basisDerFace[f], NULL);CHKERRQ(ierr); 378 } else if (id == PETSCFV_CLASSID) { 379 PetscFV fv = (PetscFV) obj; 380 381 ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr); 382 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 383 Nb = Nc; 384 ierr = PetscFVGetDefaultTabulation(fv, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr); 385 /* TODO: should PetscFV also have face tabulation? Otherwise there will be a null pointer in prob->basisFace */ 386 } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 387 prob->Nc[f] = Nc; 388 prob->Nb[f] = Nb; 389 prob->off[f+1] = Nc + prob->off[f]; 390 prob->offDer[f+1] = Nc*dim + prob->offDer[f]; 391 if (q) {ierr = PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);} 392 NqMax = PetscMax(NqMax, Nq); 393 NbMax = PetscMax(NbMax, Nb); 394 NcMax = PetscMax(NcMax, Nc); 395 prob->totDim += Nb; 396 prob->totComp += Nc; 397 /* There are two faces for all fields but the cohesive field on a hybrid cell */ 398 if (prob->isHybrid && (f < Nf-1)) prob->totDim += Nb; 399 } 400 /* Allocate works space */ 401 if (prob->isHybrid) NsMax = 2; 402 ierr = PetscMalloc3(NsMax*prob->totComp,&prob->u,NsMax*prob->totComp,&prob->u_t,NsMax*prob->totComp*dimEmbed,&prob->u_x);CHKERRQ(ierr); 403 ierr = PetscMalloc5(dimEmbed,&prob->x,NbMax*NcMax,&prob->basisReal,NbMax*NcMax*dimEmbed,&prob->basisDerReal,NbMax*NcMax,&prob->testReal,NbMax*NcMax*dimEmbed,&prob->testDerReal);CHKERRQ(ierr); 404 ierr = PetscMalloc6(NsMax*NqMax*NcMax,&prob->f0,NsMax*NqMax*NcMax*dim,&prob->f1, 405 NsMax*NsMax*NqMax*NcMax*NcMax,&prob->g0,NsMax*NsMax*NqMax*NcMax*NcMax*dim,&prob->g1, 406 NsMax*NsMax*NqMax*NcMax*NcMax*dim,&prob->g2,NsMax*NsMax*NqMax*NcMax*NcMax*dim*dim,&prob->g3);CHKERRQ(ierr); 407 if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);} 408 prob->setup = PETSC_TRUE; 409 PetscFunctionReturn(0); 410 } 411 412 static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob) 413 { 414 PetscErrorCode ierr; 415 416 PetscFunctionBegin; 417 ierr = PetscFree2(prob->Nc,prob->Nb);CHKERRQ(ierr); 418 ierr = PetscFree2(prob->off,prob->offDer);CHKERRQ(ierr); 419 ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisFace,prob->basisDerFace);CHKERRQ(ierr); 420 ierr = PetscFree3(prob->u,prob->u_t,prob->u_x);CHKERRQ(ierr); 421 ierr = PetscFree5(prob->x,prob->basisReal, prob->basisDerReal,prob->testReal,prob->testDerReal);CHKERRQ(ierr); 422 ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr); 423 PetscFunctionReturn(0); 424 } 425 426 static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew) 427 { 428 PetscObject *tmpd; 429 PetscBool *tmpi; 430 PetscPointFunc *tmpobj, *tmpf, *tmpup; 431 PetscPointJac *tmpg, *tmpgp, *tmpgt; 432 PetscBdPointFunc *tmpfbd; 433 PetscBdPointJac *tmpgbd; 434 PetscRiemannFunc *tmpr; 435 PetscSimplePointFunc *tmpexactSol; 436 void **tmpexactCtx; 437 void **tmpctx; 438 PetscInt Nf = prob->Nf, f; 439 PetscErrorCode ierr; 440 441 PetscFunctionBegin; 442 if (Nf >= NfNew) PetscFunctionReturn(0); 443 prob->setup = PETSC_FALSE; 444 ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr); 445 ierr = PetscMalloc2(NfNew, &tmpd, NfNew, &tmpi);CHKERRQ(ierr); 446 for (f = 0; f < Nf; ++f) {tmpd[f] = prob->disc[f]; tmpi[f] = prob->implicit[f];} 447 for (f = Nf; f < NfNew; ++f) {tmpd[f] = NULL; tmpi[f] = PETSC_TRUE;} 448 ierr = PetscFree2(prob->disc, prob->implicit);CHKERRQ(ierr); 449 prob->Nf = NfNew; 450 prob->disc = tmpd; 451 prob->implicit = tmpi; 452 ierr = PetscCalloc7(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew*NfNew*4, &tmpgp, NfNew*NfNew*4, &tmpgt, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr); 453 ierr = PetscCalloc1(NfNew, &tmpup);CHKERRQ(ierr); 454 for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f]; 455 for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f]; 456 for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f]; 457 for (f = 0; f < Nf*Nf*4; ++f) tmpgp[f] = prob->gp[f]; 458 for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f]; 459 for (f = 0; f < Nf; ++f) tmpup[f] = prob->update[f]; 460 for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f]; 461 for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL; 462 for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL; 463 for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL; 464 for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgp[f] = NULL; 465 for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgt[f] = NULL; 466 for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL; 467 for (f = Nf; f < NfNew; ++f) tmpup[f] = NULL; 468 for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL; 469 ierr = PetscFree7(prob->obj, prob->f, prob->g, prob->gp, prob->gt, prob->r, prob->ctx);CHKERRQ(ierr); 470 ierr = PetscFree(prob->update);CHKERRQ(ierr); 471 prob->obj = tmpobj; 472 prob->f = tmpf; 473 prob->g = tmpg; 474 prob->gp = tmpgp; 475 prob->gt = tmpgt; 476 prob->r = tmpr; 477 prob->update = tmpup; 478 prob->ctx = tmpctx; 479 ierr = PetscCalloc4(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd, NfNew, &tmpexactSol, NfNew, &tmpexactCtx);CHKERRQ(ierr); 480 for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f]; 481 for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f]; 482 for (f = 0; f < Nf; ++f) tmpexactSol[f] = prob->exactSol[f]; 483 for (f = 0; f < Nf; ++f) tmpexactCtx[f] = prob->exactCtx[f]; 484 for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL; 485 for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL; 486 for (f = Nf; f < NfNew; ++f) tmpexactSol[f] = NULL; 487 for (f = Nf; f < NfNew; ++f) tmpexactCtx[f] = NULL; 488 ierr = PetscFree4(prob->fBd, prob->gBd, prob->exactSol, prob->exactCtx);CHKERRQ(ierr); 489 prob->fBd = tmpfbd; 490 prob->gBd = tmpgbd; 491 prob->exactSol = tmpexactSol; 492 prob->exactCtx = tmpexactCtx; 493 PetscFunctionReturn(0); 494 } 495 496 /*@ 497 PetscDSDestroy - Destroys a PetscDS object 498 499 Collective on PetscDS 500 501 Input Parameter: 502 . prob - the PetscDS object to destroy 503 504 Level: developer 505 506 .seealso PetscDSView() 507 @*/ 508 PetscErrorCode PetscDSDestroy(PetscDS *prob) 509 { 510 PetscInt f; 511 DSBoundary next; 512 PetscErrorCode ierr; 513 514 PetscFunctionBegin; 515 if (!*prob) PetscFunctionReturn(0); 516 PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1); 517 518 if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);} 519 ((PetscObject) (*prob))->refct = 0; 520 if ((*prob)->subprobs) { 521 PetscInt dim, d; 522 523 ierr = PetscDSGetSpatialDimension(*prob, &dim);CHKERRQ(ierr); 524 for (d = 0; d < dim; ++d) {ierr = PetscDSDestroy(&(*prob)->subprobs[d]);CHKERRQ(ierr);} 525 } 526 ierr = PetscFree((*prob)->subprobs);CHKERRQ(ierr); 527 ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr); 528 for (f = 0; f < (*prob)->Nf; ++f) { 529 ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr); 530 } 531 ierr = PetscFree2((*prob)->disc, (*prob)->implicit);CHKERRQ(ierr); 532 ierr = PetscFree7((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->gp,(*prob)->gt,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr); 533 ierr = PetscFree((*prob)->update);CHKERRQ(ierr); 534 ierr = PetscFree4((*prob)->fBd,(*prob)->gBd,(*prob)->exactSol,(*prob)->exactCtx);CHKERRQ(ierr); 535 if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);} 536 next = (*prob)->boundary; 537 while (next) { 538 DSBoundary b = next; 539 540 next = b->next; 541 ierr = PetscFree(b->comps);CHKERRQ(ierr); 542 ierr = PetscFree(b->ids);CHKERRQ(ierr); 543 ierr = PetscFree(b->name);CHKERRQ(ierr); 544 ierr = PetscFree(b->labelname);CHKERRQ(ierr); 545 ierr = PetscFree(b);CHKERRQ(ierr); 546 } 547 ierr = PetscFree((*prob)->constants);CHKERRQ(ierr); 548 ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr); 549 PetscFunctionReturn(0); 550 } 551 552 /*@ 553 PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType(). 554 555 Collective on MPI_Comm 556 557 Input Parameter: 558 . comm - The communicator for the PetscDS object 559 560 Output Parameter: 561 . prob - The PetscDS object 562 563 Level: beginner 564 565 .seealso: PetscDSSetType(), PETSCDSBASIC 566 @*/ 567 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob) 568 { 569 PetscDS p; 570 PetscErrorCode ierr; 571 572 PetscFunctionBegin; 573 PetscValidPointer(prob, 2); 574 *prob = NULL; 575 ierr = PetscDSInitializePackage();CHKERRQ(ierr); 576 577 ierr = PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr); 578 579 p->Nf = 0; 580 p->setup = PETSC_FALSE; 581 p->numConstants = 0; 582 p->constants = NULL; 583 p->dimEmbed = -1; 584 p->useJacPre = PETSC_TRUE; 585 586 *prob = p; 587 PetscFunctionReturn(0); 588 } 589 590 /*@ 591 PetscDSGetNumFields - Returns the number of fields in the DS 592 593 Not collective 594 595 Input Parameter: 596 . prob - The PetscDS object 597 598 Output Parameter: 599 . Nf - The number of fields 600 601 Level: beginner 602 603 .seealso: PetscDSGetSpatialDimension(), PetscDSCreate() 604 @*/ 605 PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf) 606 { 607 PetscFunctionBegin; 608 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 609 PetscValidPointer(Nf, 2); 610 *Nf = prob->Nf; 611 PetscFunctionReturn(0); 612 } 613 614 /*@ 615 PetscDSGetSpatialDimension - Returns the spatial dimension of the DS, meaning the topological dimension of the discretizations 616 617 Not collective 618 619 Input Parameter: 620 . prob - The PetscDS object 621 622 Output Parameter: 623 . dim - The spatial dimension 624 625 Level: beginner 626 627 .seealso: PetscDSGetCoordinateDimension(), PetscDSGetNumFields(), PetscDSCreate() 628 @*/ 629 PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim) 630 { 631 PetscErrorCode ierr; 632 633 PetscFunctionBegin; 634 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 635 PetscValidPointer(dim, 2); 636 *dim = 0; 637 if (prob->Nf) { 638 PetscObject obj; 639 PetscClassId id; 640 641 ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr); 642 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 643 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);} 644 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);} 645 else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0); 646 } 647 PetscFunctionReturn(0); 648 } 649 650 /*@ 651 PetscDSGetCoordinateDimension - Returns the coordinate dimension of the DS, meaning the dimension of the space into which the discretiaztions are embedded 652 653 Not collective 654 655 Input Parameter: 656 . prob - The PetscDS object 657 658 Output Parameter: 659 . dimEmbed - The coordinate dimension 660 661 Level: beginner 662 663 .seealso: PetscDSSetCoordinateDimension(), PetscDSGetSpatialDimension(), PetscDSGetNumFields(), PetscDSCreate() 664 @*/ 665 PetscErrorCode PetscDSGetCoordinateDimension(PetscDS prob, PetscInt *dimEmbed) 666 { 667 PetscFunctionBegin; 668 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 669 PetscValidPointer(dimEmbed, 2); 670 if (prob->dimEmbed < 0) SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONGSTATE, "No coordinate dimension set for this DS"); 671 *dimEmbed = prob->dimEmbed; 672 PetscFunctionReturn(0); 673 } 674 675 /*@ 676 PetscDSSetCoordinateDimension - Set the coordinate dimension of the DS, meaning the dimension of the space into which the discretiaztions are embedded 677 678 Logically collective on DS 679 680 Input Parameters: 681 + prob - The PetscDS object 682 - dimEmbed - The coordinate dimension 683 684 Level: beginner 685 686 .seealso: PetscDSGetCoordinateDimension(), PetscDSGetSpatialDimension(), PetscDSGetNumFields(), PetscDSCreate() 687 @*/ 688 PetscErrorCode PetscDSSetCoordinateDimension(PetscDS prob, PetscInt dimEmbed) 689 { 690 PetscFunctionBegin; 691 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 692 if (dimEmbed < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coordinate dimension must be non-negative, not %D", dimEmbed); 693 prob->dimEmbed = dimEmbed; 694 PetscFunctionReturn(0); 695 } 696 697 /*@ 698 PetscDSGetHybrid - Returns the flag for a hybrid (cohesive) cell 699 700 Not collective 701 702 Input Parameter: 703 . prob - The PetscDS object 704 705 Output Parameter: 706 . isHybrid - The flag 707 708 Level: developer 709 710 .seealso: PetscDSSetHybrid(), PetscDSCreate() 711 @*/ 712 PetscErrorCode PetscDSGetHybrid(PetscDS prob, PetscBool *isHybrid) 713 { 714 PetscFunctionBegin; 715 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 716 PetscValidPointer(isHybrid, 2); 717 *isHybrid = prob->isHybrid; 718 PetscFunctionReturn(0); 719 } 720 721 /*@ 722 PetscDSSetHybrid - Set the flag for a hybrid (cohesive) cell 723 724 Not collective 725 726 Input Parameters: 727 + prob - The PetscDS object 728 - isHybrid - The flag 729 730 Level: developer 731 732 .seealso: PetscDSGetHybrid(), PetscDSCreate() 733 @*/ 734 PetscErrorCode PetscDSSetHybrid(PetscDS prob, PetscBool isHybrid) 735 { 736 PetscFunctionBegin; 737 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 738 prob->isHybrid = isHybrid; 739 PetscFunctionReturn(0); 740 } 741 742 /*@ 743 PetscDSGetTotalDimension - Returns the total size of the approximation space for this system 744 745 Not collective 746 747 Input Parameter: 748 . prob - The PetscDS object 749 750 Output Parameter: 751 . dim - The total problem dimension 752 753 Level: beginner 754 755 .seealso: PetscDSGetNumFields(), PetscDSCreate() 756 @*/ 757 PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim) 758 { 759 PetscErrorCode ierr; 760 761 PetscFunctionBegin; 762 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 763 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 764 PetscValidPointer(dim, 2); 765 *dim = prob->totDim; 766 PetscFunctionReturn(0); 767 } 768 769 /*@ 770 PetscDSGetTotalComponents - Returns the total number of components in this system 771 772 Not collective 773 774 Input Parameter: 775 . prob - The PetscDS object 776 777 Output Parameter: 778 . dim - The total number of components 779 780 Level: beginner 781 782 .seealso: PetscDSGetNumFields(), PetscDSCreate() 783 @*/ 784 PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc) 785 { 786 PetscErrorCode ierr; 787 788 PetscFunctionBegin; 789 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 790 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 791 PetscValidPointer(Nc, 2); 792 *Nc = prob->totComp; 793 PetscFunctionReturn(0); 794 } 795 796 /*@ 797 PetscDSGetDiscretization - Returns the discretization object for the given field 798 799 Not collective 800 801 Input Parameters: 802 + prob - The PetscDS object 803 - f - The field number 804 805 Output Parameter: 806 . disc - The discretization object 807 808 Level: beginner 809 810 .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 811 @*/ 812 PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc) 813 { 814 PetscFunctionBegin; 815 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 816 PetscValidPointer(disc, 3); 817 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 818 *disc = prob->disc[f]; 819 PetscFunctionReturn(0); 820 } 821 822 /*@ 823 PetscDSSetDiscretization - Sets the discretization object for the given field 824 825 Not collective 826 827 Input Parameters: 828 + prob - The PetscDS object 829 . f - The field number 830 - disc - The discretization object 831 832 Level: beginner 833 834 .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 835 @*/ 836 PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc) 837 { 838 PetscErrorCode ierr; 839 840 PetscFunctionBegin; 841 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 842 PetscValidPointer(disc, 3); 843 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 844 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 845 ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr); 846 prob->disc[f] = disc; 847 ierr = PetscObjectReference(disc);CHKERRQ(ierr); 848 { 849 PetscClassId id; 850 851 ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr); 852 if (id == PETSCFE_CLASSID) { 853 ierr = PetscDSSetImplicit(prob, f, PETSC_TRUE);CHKERRQ(ierr); 854 } else if (id == PETSCFV_CLASSID) { 855 ierr = PetscDSSetImplicit(prob, f, PETSC_FALSE);CHKERRQ(ierr); 856 } 857 } 858 PetscFunctionReturn(0); 859 } 860 861 /*@ 862 PetscDSAddDiscretization - Adds a discretization object 863 864 Not collective 865 866 Input Parameters: 867 + prob - The PetscDS object 868 - disc - The boundary discretization object 869 870 Level: beginner 871 872 .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 873 @*/ 874 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc) 875 { 876 PetscErrorCode ierr; 877 878 PetscFunctionBegin; 879 ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr); 880 PetscFunctionReturn(0); 881 } 882 883 /*@ 884 PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX 885 886 Not collective 887 888 Input Parameters: 889 + prob - The PetscDS object 890 - f - The field number 891 892 Output Parameter: 893 . implicit - The flag indicating what kind of solve to use for this field 894 895 Level: developer 896 897 .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 898 @*/ 899 PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit) 900 { 901 PetscFunctionBegin; 902 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 903 PetscValidPointer(implicit, 3); 904 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 905 *implicit = prob->implicit[f]; 906 PetscFunctionReturn(0); 907 } 908 909 /*@ 910 PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX 911 912 Not collective 913 914 Input Parameters: 915 + prob - The PetscDS object 916 . f - The field number 917 - implicit - The flag indicating what kind of solve to use for this field 918 919 Level: developer 920 921 .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 922 @*/ 923 PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit) 924 { 925 PetscFunctionBegin; 926 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 927 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 928 prob->implicit[f] = implicit; 929 PetscFunctionReturn(0); 930 } 931 932 PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f, 933 void (**obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 934 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 935 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 936 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])) 937 { 938 PetscFunctionBegin; 939 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 940 PetscValidPointer(obj, 2); 941 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 942 *obj = prob->obj[f]; 943 PetscFunctionReturn(0); 944 } 945 946 PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f, 947 void (*obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 948 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 949 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 950 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])) 951 { 952 PetscErrorCode ierr; 953 954 PetscFunctionBegin; 955 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 956 if (obj) PetscValidFunction(obj, 2); 957 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 958 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 959 prob->obj[f] = obj; 960 PetscFunctionReturn(0); 961 } 962 963 /*@C 964 PetscDSGetResidual - Get the pointwise residual function for a given test field 965 966 Not collective 967 968 Input Parameters: 969 + prob - The PetscDS 970 - f - The test field number 971 972 Output Parameters: 973 + f0 - integrand for the test function term 974 - f1 - integrand for the test function gradient term 975 976 Note: We are using a first order FEM model for the weak form: 977 978 \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) 979 980 The calling sequence for the callbacks f0 and f1 is given by: 981 982 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 983 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 984 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 985 $ PetscReal t, const PetscReal x[], PetscScalar f0[]) 986 987 + dim - the spatial dimension 988 . Nf - the number of fields 989 . uOff - the offset into u[] and u_t[] for each field 990 . uOff_x - the offset into u_x[] for each field 991 . u - each field evaluated at the current point 992 . u_t - the time derivative of each field evaluated at the current point 993 . u_x - the gradient of each field evaluated at the current point 994 . aOff - the offset into a[] and a_t[] for each auxiliary field 995 . aOff_x - the offset into a_x[] for each auxiliary field 996 . a - each auxiliary field evaluated at the current point 997 . a_t - the time derivative of each auxiliary field evaluated at the current point 998 . a_x - the gradient of auxiliary each field evaluated at the current point 999 . t - current time 1000 . x - coordinates of the current point 1001 . numConstants - number of constant parameters 1002 . constants - constant parameters 1003 - f0 - output values at the current point 1004 1005 Level: intermediate 1006 1007 .seealso: PetscDSSetResidual() 1008 @*/ 1009 PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f, 1010 void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1011 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1012 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1013 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), 1014 void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1015 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1016 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1017 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 1018 { 1019 PetscFunctionBegin; 1020 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1021 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1022 if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];} 1023 if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];} 1024 PetscFunctionReturn(0); 1025 } 1026 1027 /*@C 1028 PetscDSSetResidual - Set the pointwise residual function for a given test field 1029 1030 Not collective 1031 1032 Input Parameters: 1033 + prob - The PetscDS 1034 . f - The test field number 1035 . f0 - integrand for the test function term 1036 - f1 - integrand for the test function gradient term 1037 1038 Note: We are using a first order FEM model for the weak form: 1039 1040 \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) 1041 1042 The calling sequence for the callbacks f0 and f1 is given by: 1043 1044 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1045 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1046 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1047 $ PetscReal t, const PetscReal x[], PetscScalar f0[]) 1048 1049 + dim - the spatial dimension 1050 . Nf - the number of fields 1051 . uOff - the offset into u[] and u_t[] for each field 1052 . uOff_x - the offset into u_x[] for each field 1053 . u - each field evaluated at the current point 1054 . u_t - the time derivative of each field evaluated at the current point 1055 . u_x - the gradient of each field evaluated at the current point 1056 . aOff - the offset into a[] and a_t[] for each auxiliary field 1057 . aOff_x - the offset into a_x[] for each auxiliary field 1058 . a - each auxiliary field evaluated at the current point 1059 . a_t - the time derivative of each auxiliary field evaluated at the current point 1060 . a_x - the gradient of auxiliary each field evaluated at the current point 1061 . t - current time 1062 . x - coordinates of the current point 1063 . numConstants - number of constant parameters 1064 . constants - constant parameters 1065 - f0 - output values at the current point 1066 1067 Level: intermediate 1068 1069 .seealso: PetscDSGetResidual() 1070 @*/ 1071 PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f, 1072 void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1073 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1074 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1075 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), 1076 void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1077 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1078 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1079 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 1080 { 1081 PetscErrorCode ierr; 1082 1083 PetscFunctionBegin; 1084 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1085 if (f0) PetscValidFunction(f0, 3); 1086 if (f1) PetscValidFunction(f1, 4); 1087 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1088 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 1089 prob->f[f*2+0] = f0; 1090 prob->f[f*2+1] = f1; 1091 PetscFunctionReturn(0); 1092 } 1093 1094 /*@C 1095 PetscDSHasJacobian - Signals that Jacobian functions have been set 1096 1097 Not collective 1098 1099 Input Parameter: 1100 . prob - The PetscDS 1101 1102 Output Parameter: 1103 . hasJac - flag that pointwise function for the Jacobian has been set 1104 1105 Level: intermediate 1106 1107 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian() 1108 @*/ 1109 PetscErrorCode PetscDSHasJacobian(PetscDS prob, PetscBool *hasJac) 1110 { 1111 PetscInt f, g, h; 1112 1113 PetscFunctionBegin; 1114 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1115 *hasJac = PETSC_FALSE; 1116 for (f = 0; f < prob->Nf; ++f) { 1117 for (g = 0; g < prob->Nf; ++g) { 1118 for (h = 0; h < 4; ++h) { 1119 if (prob->g[(f*prob->Nf + g)*4+h]) *hasJac = PETSC_TRUE; 1120 } 1121 } 1122 } 1123 PetscFunctionReturn(0); 1124 } 1125 1126 /*@C 1127 PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field 1128 1129 Not collective 1130 1131 Input Parameters: 1132 + prob - The PetscDS 1133 . f - The test field number 1134 - g - The field number 1135 1136 Output Parameters: 1137 + g0 - integrand for the test and basis function term 1138 . g1 - integrand for the test function and basis function gradient term 1139 . g2 - integrand for the test function gradient and basis function term 1140 - g3 - integrand for the test function gradient and basis function gradient term 1141 1142 Note: We are using a first order FEM model for the weak form: 1143 1144 \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 1145 1146 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1147 1148 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1149 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1150 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1151 $ PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]) 1152 1153 + dim - the spatial dimension 1154 . Nf - the number of fields 1155 . uOff - the offset into u[] and u_t[] for each field 1156 . uOff_x - the offset into u_x[] for each field 1157 . u - each field evaluated at the current point 1158 . u_t - the time derivative of each field evaluated at the current point 1159 . u_x - the gradient of each field evaluated at the current point 1160 . aOff - the offset into a[] and a_t[] for each auxiliary field 1161 . aOff_x - the offset into a_x[] for each auxiliary field 1162 . a - each auxiliary field evaluated at the current point 1163 . a_t - the time derivative of each auxiliary field evaluated at the current point 1164 . a_x - the gradient of auxiliary each field evaluated at the current point 1165 . t - current time 1166 . u_tShift - the multiplier a for dF/dU_t 1167 . x - coordinates of the current point 1168 . numConstants - number of constant parameters 1169 . constants - constant parameters 1170 - g0 - output values at the current point 1171 1172 Level: intermediate 1173 1174 .seealso: PetscDSSetJacobian() 1175 @*/ 1176 PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g, 1177 void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1178 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1179 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1180 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1181 void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1182 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1183 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1184 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1185 void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1186 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1187 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1188 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1189 void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1190 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1191 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1192 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1193 { 1194 PetscFunctionBegin; 1195 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1196 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1197 if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf); 1198 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];} 1199 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];} 1200 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];} 1201 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];} 1202 PetscFunctionReturn(0); 1203 } 1204 1205 /*@C 1206 PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields 1207 1208 Not collective 1209 1210 Input Parameters: 1211 + prob - The PetscDS 1212 . f - The test field number 1213 . g - The field number 1214 . g0 - integrand for the test and basis function term 1215 . g1 - integrand for the test function and basis function gradient term 1216 . g2 - integrand for the test function gradient and basis function term 1217 - g3 - integrand for the test function gradient and basis function gradient term 1218 1219 Note: We are using a first order FEM model for the weak form: 1220 1221 \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 1222 1223 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1224 1225 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1226 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1227 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1228 $ PetscReal t, const PetscReal x[], PetscScalar g0[]) 1229 1230 + dim - the spatial dimension 1231 . Nf - the number of fields 1232 . uOff - the offset into u[] and u_t[] for each field 1233 . uOff_x - the offset into u_x[] for each field 1234 . u - each field evaluated at the current point 1235 . u_t - the time derivative of each field evaluated at the current point 1236 . u_x - the gradient of each field evaluated at the current point 1237 . aOff - the offset into a[] and a_t[] for each auxiliary field 1238 . aOff_x - the offset into a_x[] for each auxiliary field 1239 . a - each auxiliary field evaluated at the current point 1240 . a_t - the time derivative of each auxiliary field evaluated at the current point 1241 . a_x - the gradient of auxiliary each field evaluated at the current point 1242 . t - current time 1243 . u_tShift - the multiplier a for dF/dU_t 1244 . x - coordinates of the current point 1245 . numConstants - number of constant parameters 1246 . constants - constant parameters 1247 - g0 - output values at the current point 1248 1249 Level: intermediate 1250 1251 .seealso: PetscDSGetJacobian() 1252 @*/ 1253 PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g, 1254 void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1255 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1256 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1257 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1258 void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1259 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1260 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1261 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1262 void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1263 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1264 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1265 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1266 void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1267 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1268 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1269 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1270 { 1271 PetscErrorCode ierr; 1272 1273 PetscFunctionBegin; 1274 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1275 if (g0) PetscValidFunction(g0, 4); 1276 if (g1) PetscValidFunction(g1, 5); 1277 if (g2) PetscValidFunction(g2, 6); 1278 if (g3) PetscValidFunction(g3, 7); 1279 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1280 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 1281 ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr); 1282 prob->g[(f*prob->Nf + g)*4+0] = g0; 1283 prob->g[(f*prob->Nf + g)*4+1] = g1; 1284 prob->g[(f*prob->Nf + g)*4+2] = g2; 1285 prob->g[(f*prob->Nf + g)*4+3] = g3; 1286 PetscFunctionReturn(0); 1287 } 1288 1289 /*@C 1290 PetscDSUseJacobianPreconditioner - Whether to construct a Jacobian preconditioner 1291 1292 Not collective 1293 1294 Input Parameters: 1295 + prob - The PetscDS 1296 - useJacPre - flag that enables construction of a Jacobian preconditioner 1297 1298 Level: intermediate 1299 1300 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian() 1301 @*/ 1302 PetscErrorCode PetscDSUseJacobianPreconditioner(PetscDS prob, PetscBool useJacPre) 1303 { 1304 PetscFunctionBegin; 1305 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1306 prob->useJacPre = useJacPre; 1307 PetscFunctionReturn(0); 1308 } 1309 1310 /*@C 1311 PetscDSHasJacobianPreconditioner - Signals that a Jacobian preconditioner matrix has been set 1312 1313 Not collective 1314 1315 Input Parameter: 1316 . prob - The PetscDS 1317 1318 Output Parameter: 1319 . hasJacPre - flag that pointwise function for Jacobian preconditioner matrix has been set 1320 1321 Level: intermediate 1322 1323 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian() 1324 @*/ 1325 PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS prob, PetscBool *hasJacPre) 1326 { 1327 PetscInt f, g, h; 1328 1329 PetscFunctionBegin; 1330 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1331 *hasJacPre = PETSC_FALSE; 1332 if (!prob->useJacPre) PetscFunctionReturn(0); 1333 for (f = 0; f < prob->Nf; ++f) { 1334 for (g = 0; g < prob->Nf; ++g) { 1335 for (h = 0; h < 4; ++h) { 1336 if (prob->gp[(f*prob->Nf + g)*4+h]) *hasJacPre = PETSC_TRUE; 1337 } 1338 } 1339 } 1340 PetscFunctionReturn(0); 1341 } 1342 1343 /*@C 1344 PetscDSGetJacobianPreconditioner - Get the pointwise Jacobian preconditioner function for given test and basis field. If this is missing, the system matrix is used to build the preconditioner. 1345 1346 Not collective 1347 1348 Input Parameters: 1349 + prob - The PetscDS 1350 . f - The test field number 1351 - g - The field number 1352 1353 Output Parameters: 1354 + g0 - integrand for the test and basis function term 1355 . g1 - integrand for the test function and basis function gradient term 1356 . g2 - integrand for the test function gradient and basis function term 1357 - g3 - integrand for the test function gradient and basis function gradient term 1358 1359 Note: We are using a first order FEM model for the weak form: 1360 1361 \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 1362 1363 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1364 1365 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1366 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1367 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1368 $ PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]) 1369 1370 + dim - the spatial dimension 1371 . Nf - the number of fields 1372 . uOff - the offset into u[] and u_t[] for each field 1373 . uOff_x - the offset into u_x[] for each field 1374 . u - each field evaluated at the current point 1375 . u_t - the time derivative of each field evaluated at the current point 1376 . u_x - the gradient of each field evaluated at the current point 1377 . aOff - the offset into a[] and a_t[] for each auxiliary field 1378 . aOff_x - the offset into a_x[] for each auxiliary field 1379 . a - each auxiliary field evaluated at the current point 1380 . a_t - the time derivative of each auxiliary field evaluated at the current point 1381 . a_x - the gradient of auxiliary each field evaluated at the current point 1382 . t - current time 1383 . u_tShift - the multiplier a for dF/dU_t 1384 . x - coordinates of the current point 1385 . numConstants - number of constant parameters 1386 . constants - constant parameters 1387 - g0 - output values at the current point 1388 1389 Level: intermediate 1390 1391 .seealso: PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian() 1392 @*/ 1393 PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS prob, PetscInt f, PetscInt g, 1394 void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1395 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1396 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1397 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1398 void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1399 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1400 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1401 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1402 void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1403 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1404 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1405 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1406 void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1407 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1408 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1409 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1410 { 1411 PetscFunctionBegin; 1412 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1413 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1414 if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf); 1415 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gp[(f*prob->Nf + g)*4+0];} 1416 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gp[(f*prob->Nf + g)*4+1];} 1417 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gp[(f*prob->Nf + g)*4+2];} 1418 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gp[(f*prob->Nf + g)*4+3];} 1419 PetscFunctionReturn(0); 1420 } 1421 1422 /*@C 1423 PetscDSSetJacobianPreconditioner - Set the pointwise Jacobian preconditioner function for given test and basis fields. If this is missing, the system matrix is used to build the preconditioner. 1424 1425 Not collective 1426 1427 Input Parameters: 1428 + prob - The PetscDS 1429 . f - The test field number 1430 . g - The field number 1431 . g0 - integrand for the test and basis function term 1432 . g1 - integrand for the test function and basis function gradient term 1433 . g2 - integrand for the test function gradient and basis function term 1434 - g3 - integrand for the test function gradient and basis function gradient term 1435 1436 Note: We are using a first order FEM model for the weak form: 1437 1438 \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 1439 1440 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1441 1442 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1443 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1444 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1445 $ PetscReal t, const PetscReal x[], PetscScalar g0[]) 1446 1447 + dim - the spatial dimension 1448 . Nf - the number of fields 1449 . uOff - the offset into u[] and u_t[] for each field 1450 . uOff_x - the offset into u_x[] for each field 1451 . u - each field evaluated at the current point 1452 . u_t - the time derivative of each field evaluated at the current point 1453 . u_x - the gradient of each field evaluated at the current point 1454 . aOff - the offset into a[] and a_t[] for each auxiliary field 1455 . aOff_x - the offset into a_x[] for each auxiliary field 1456 . a - each auxiliary field evaluated at the current point 1457 . a_t - the time derivative of each auxiliary field evaluated at the current point 1458 . a_x - the gradient of auxiliary each field evaluated at the current point 1459 . t - current time 1460 . u_tShift - the multiplier a for dF/dU_t 1461 . x - coordinates of the current point 1462 . numConstants - number of constant parameters 1463 . constants - constant parameters 1464 - g0 - output values at the current point 1465 1466 Level: intermediate 1467 1468 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobian() 1469 @*/ 1470 PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS prob, PetscInt f, PetscInt g, 1471 void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1472 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1473 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1474 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1475 void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1476 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1477 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1478 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1479 void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1480 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1481 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1482 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1483 void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1484 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1485 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1486 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1487 { 1488 PetscErrorCode ierr; 1489 1490 PetscFunctionBegin; 1491 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1492 if (g0) PetscValidFunction(g0, 4); 1493 if (g1) PetscValidFunction(g1, 5); 1494 if (g2) PetscValidFunction(g2, 6); 1495 if (g3) PetscValidFunction(g3, 7); 1496 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1497 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 1498 ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr); 1499 prob->gp[(f*prob->Nf + g)*4+0] = g0; 1500 prob->gp[(f*prob->Nf + g)*4+1] = g1; 1501 prob->gp[(f*prob->Nf + g)*4+2] = g2; 1502 prob->gp[(f*prob->Nf + g)*4+3] = g3; 1503 PetscFunctionReturn(0); 1504 } 1505 1506 /*@C 1507 PetscDSHasDynamicJacobian - Signals that a dynamic Jacobian, dF/du_t, has been set 1508 1509 Not collective 1510 1511 Input Parameter: 1512 . prob - The PetscDS 1513 1514 Output Parameter: 1515 . hasDynJac - flag that pointwise function for dynamic Jacobian has been set 1516 1517 Level: intermediate 1518 1519 .seealso: PetscDSGetDynamicJacobian(), PetscDSSetDynamicJacobian(), PetscDSGetJacobian() 1520 @*/ 1521 PetscErrorCode PetscDSHasDynamicJacobian(PetscDS prob, PetscBool *hasDynJac) 1522 { 1523 PetscInt f, g, h; 1524 1525 PetscFunctionBegin; 1526 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1527 *hasDynJac = PETSC_FALSE; 1528 for (f = 0; f < prob->Nf; ++f) { 1529 for (g = 0; g < prob->Nf; ++g) { 1530 for (h = 0; h < 4; ++h) { 1531 if (prob->gt[(f*prob->Nf + g)*4+h]) *hasDynJac = PETSC_TRUE; 1532 } 1533 } 1534 } 1535 PetscFunctionReturn(0); 1536 } 1537 1538 /*@C 1539 PetscDSGetDynamicJacobian - Get the pointwise dynamic Jacobian, dF/du_t, function for given test and basis field 1540 1541 Not collective 1542 1543 Input Parameters: 1544 + prob - The PetscDS 1545 . f - The test field number 1546 - g - The field number 1547 1548 Output Parameters: 1549 + g0 - integrand for the test and basis function term 1550 . g1 - integrand for the test function and basis function gradient term 1551 . g2 - integrand for the test function gradient and basis function term 1552 - g3 - integrand for the test function gradient and basis function gradient term 1553 1554 Note: We are using a first order FEM model for the weak form: 1555 1556 \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 1557 1558 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1559 1560 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1561 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1562 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1563 $ PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]) 1564 1565 + dim - the spatial dimension 1566 . Nf - the number of fields 1567 . uOff - the offset into u[] and u_t[] for each field 1568 . uOff_x - the offset into u_x[] for each field 1569 . u - each field evaluated at the current point 1570 . u_t - the time derivative of each field evaluated at the current point 1571 . u_x - the gradient of each field evaluated at the current point 1572 . aOff - the offset into a[] and a_t[] for each auxiliary field 1573 . aOff_x - the offset into a_x[] for each auxiliary field 1574 . a - each auxiliary field evaluated at the current point 1575 . a_t - the time derivative of each auxiliary field evaluated at the current point 1576 . a_x - the gradient of auxiliary each field evaluated at the current point 1577 . t - current time 1578 . u_tShift - the multiplier a for dF/dU_t 1579 . x - coordinates of the current point 1580 . numConstants - number of constant parameters 1581 . constants - constant parameters 1582 - g0 - output values at the current point 1583 1584 Level: intermediate 1585 1586 .seealso: PetscDSSetJacobian() 1587 @*/ 1588 PetscErrorCode PetscDSGetDynamicJacobian(PetscDS prob, PetscInt f, PetscInt g, 1589 void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1590 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1591 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1592 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1593 void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1594 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1595 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1596 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1597 void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1598 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1599 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1600 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1601 void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1602 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1603 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1604 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1605 { 1606 PetscFunctionBegin; 1607 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1608 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1609 if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf); 1610 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gt[(f*prob->Nf + g)*4+0];} 1611 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gt[(f*prob->Nf + g)*4+1];} 1612 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gt[(f*prob->Nf + g)*4+2];} 1613 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gt[(f*prob->Nf + g)*4+3];} 1614 PetscFunctionReturn(0); 1615 } 1616 1617 /*@C 1618 PetscDSSetDynamicJacobian - Set the pointwise dynamic Jacobian, dF/du_t, function for given test and basis fields 1619 1620 Not collective 1621 1622 Input Parameters: 1623 + prob - The PetscDS 1624 . f - The test field number 1625 . g - The field number 1626 . g0 - integrand for the test and basis function term 1627 . g1 - integrand for the test function and basis function gradient term 1628 . g2 - integrand for the test function gradient and basis function term 1629 - g3 - integrand for the test function gradient and basis function gradient term 1630 1631 Note: We are using a first order FEM model for the weak form: 1632 1633 \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 1634 1635 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1636 1637 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1638 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1639 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1640 $ PetscReal t, const PetscReal x[], PetscScalar g0[]) 1641 1642 + dim - the spatial dimension 1643 . Nf - the number of fields 1644 . uOff - the offset into u[] and u_t[] for each field 1645 . uOff_x - the offset into u_x[] for each field 1646 . u - each field evaluated at the current point 1647 . u_t - the time derivative of each field evaluated at the current point 1648 . u_x - the gradient of each field evaluated at the current point 1649 . aOff - the offset into a[] and a_t[] for each auxiliary field 1650 . aOff_x - the offset into a_x[] for each auxiliary field 1651 . a - each auxiliary field evaluated at the current point 1652 . a_t - the time derivative of each auxiliary field evaluated at the current point 1653 . a_x - the gradient of auxiliary each field evaluated at the current point 1654 . t - current time 1655 . u_tShift - the multiplier a for dF/dU_t 1656 . x - coordinates of the current point 1657 . numConstants - number of constant parameters 1658 . constants - constant parameters 1659 - g0 - output values at the current point 1660 1661 Level: intermediate 1662 1663 .seealso: PetscDSGetJacobian() 1664 @*/ 1665 PetscErrorCode PetscDSSetDynamicJacobian(PetscDS prob, PetscInt f, PetscInt g, 1666 void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1667 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1668 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1669 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1670 void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1671 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1672 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1673 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1674 void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1675 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1676 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1677 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1678 void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1679 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1680 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1681 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1682 { 1683 PetscErrorCode ierr; 1684 1685 PetscFunctionBegin; 1686 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1687 if (g0) PetscValidFunction(g0, 4); 1688 if (g1) PetscValidFunction(g1, 5); 1689 if (g2) PetscValidFunction(g2, 6); 1690 if (g3) PetscValidFunction(g3, 7); 1691 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1692 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 1693 ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr); 1694 prob->gt[(f*prob->Nf + g)*4+0] = g0; 1695 prob->gt[(f*prob->Nf + g)*4+1] = g1; 1696 prob->gt[(f*prob->Nf + g)*4+2] = g2; 1697 prob->gt[(f*prob->Nf + g)*4+3] = g3; 1698 PetscFunctionReturn(0); 1699 } 1700 1701 /*@C 1702 PetscDSGetRiemannSolver - Returns the Riemann solver for the given field 1703 1704 Not collective 1705 1706 Input Arguments: 1707 + prob - The PetscDS object 1708 - f - The field number 1709 1710 Output Argument: 1711 . r - Riemann solver 1712 1713 Calling sequence for r: 1714 1715 $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 1716 1717 + dim - The spatial dimension 1718 . Nf - The number of fields 1719 . x - The coordinates at a point on the interface 1720 . n - The normal vector to the interface 1721 . uL - The state vector to the left of the interface 1722 . uR - The state vector to the right of the interface 1723 . flux - output array of flux through the interface 1724 . numConstants - number of constant parameters 1725 . constants - constant parameters 1726 - ctx - optional user context 1727 1728 Level: intermediate 1729 1730 .seealso: PetscDSSetRiemannSolver() 1731 @*/ 1732 PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f, 1733 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)) 1734 { 1735 PetscFunctionBegin; 1736 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1737 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1738 PetscValidPointer(r, 3); 1739 *r = prob->r[f]; 1740 PetscFunctionReturn(0); 1741 } 1742 1743 /*@C 1744 PetscDSSetRiemannSolver - Sets the Riemann solver for the given field 1745 1746 Not collective 1747 1748 Input Arguments: 1749 + prob - The PetscDS object 1750 . f - The field number 1751 - r - Riemann solver 1752 1753 Calling sequence for r: 1754 1755 $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 1756 1757 + dim - The spatial dimension 1758 . Nf - The number of fields 1759 . x - The coordinates at a point on the interface 1760 . n - The normal vector to the interface 1761 . uL - The state vector to the left of the interface 1762 . uR - The state vector to the right of the interface 1763 . flux - output array of flux through the interface 1764 . numConstants - number of constant parameters 1765 . constants - constant parameters 1766 - ctx - optional user context 1767 1768 Level: intermediate 1769 1770 .seealso: PetscDSGetRiemannSolver() 1771 @*/ 1772 PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f, 1773 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)) 1774 { 1775 PetscErrorCode ierr; 1776 1777 PetscFunctionBegin; 1778 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1779 if (r) PetscValidFunction(r, 3); 1780 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1781 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 1782 prob->r[f] = r; 1783 PetscFunctionReturn(0); 1784 } 1785 1786 /*@C 1787 PetscDSGetUpdate - Get the pointwise update function for a given field 1788 1789 Not collective 1790 1791 Input Parameters: 1792 + prob - The PetscDS 1793 - f - The field number 1794 1795 Output Parameters: 1796 + update - update function 1797 1798 Note: The calling sequence for the callback update is given by: 1799 1800 $ update(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1801 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1802 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1803 $ PetscReal t, const PetscReal x[], PetscScalar uNew[]) 1804 1805 + dim - the spatial dimension 1806 . Nf - the number of fields 1807 . uOff - the offset into u[] and u_t[] for each field 1808 . uOff_x - the offset into u_x[] for each field 1809 . u - each field evaluated at the current point 1810 . u_t - the time derivative of each field evaluated at the current point 1811 . u_x - the gradient of each field evaluated at the current point 1812 . aOff - the offset into a[] and a_t[] for each auxiliary field 1813 . aOff_x - the offset into a_x[] for each auxiliary field 1814 . a - each auxiliary field evaluated at the current point 1815 . a_t - the time derivative of each auxiliary field evaluated at the current point 1816 . a_x - the gradient of auxiliary each field evaluated at the current point 1817 . t - current time 1818 . x - coordinates of the current point 1819 - uNew - new value for field at the current point 1820 1821 Level: intermediate 1822 1823 .seealso: PetscDSSetUpdate(), PetscDSSetResidual() 1824 @*/ 1825 PetscErrorCode PetscDSGetUpdate(PetscDS prob, PetscInt f, 1826 void (**update)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1827 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1828 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1829 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uNew[])) 1830 { 1831 PetscFunctionBegin; 1832 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1833 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1834 if (update) {PetscValidPointer(update, 3); *update = prob->update[f];} 1835 PetscFunctionReturn(0); 1836 } 1837 1838 /*@C 1839 PetscDSSetUpdate - Set the pointwise update function for a given field 1840 1841 Not collective 1842 1843 Input Parameters: 1844 + prob - The PetscDS 1845 . f - The field number 1846 - update - update function 1847 1848 Note: The calling sequence for the callback update is given by: 1849 1850 $ update(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1851 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1852 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1853 $ PetscReal t, const PetscReal x[], PetscScalar uNew[]) 1854 1855 + dim - the spatial dimension 1856 . Nf - the number of fields 1857 . uOff - the offset into u[] and u_t[] for each field 1858 . uOff_x - the offset into u_x[] for each field 1859 . u - each field evaluated at the current point 1860 . u_t - the time derivative of each field evaluated at the current point 1861 . u_x - the gradient of each field evaluated at the current point 1862 . aOff - the offset into a[] and a_t[] for each auxiliary field 1863 . aOff_x - the offset into a_x[] for each auxiliary field 1864 . a - each auxiliary field evaluated at the current point 1865 . a_t - the time derivative of each auxiliary field evaluated at the current point 1866 . a_x - the gradient of auxiliary each field evaluated at the current point 1867 . t - current time 1868 . x - coordinates of the current point 1869 - uNew - new field values at the current point 1870 1871 Level: intermediate 1872 1873 .seealso: PetscDSGetResidual() 1874 @*/ 1875 PetscErrorCode PetscDSSetUpdate(PetscDS prob, PetscInt f, 1876 void (*update)(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 x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uNew[])) 1880 { 1881 PetscErrorCode ierr; 1882 1883 PetscFunctionBegin; 1884 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1885 if (update) PetscValidFunction(update, 3); 1886 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1887 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 1888 prob->update[f] = update; 1889 PetscFunctionReturn(0); 1890 } 1891 1892 PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx) 1893 { 1894 PetscFunctionBegin; 1895 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1896 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1897 PetscValidPointer(ctx, 3); 1898 *ctx = prob->ctx[f]; 1899 PetscFunctionReturn(0); 1900 } 1901 1902 PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx) 1903 { 1904 PetscErrorCode ierr; 1905 1906 PetscFunctionBegin; 1907 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1908 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1909 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 1910 prob->ctx[f] = ctx; 1911 PetscFunctionReturn(0); 1912 } 1913 1914 /*@C 1915 PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field 1916 1917 Not collective 1918 1919 Input Parameters: 1920 + prob - The PetscDS 1921 - f - The test field number 1922 1923 Output Parameters: 1924 + f0 - boundary integrand for the test function term 1925 - f1 - boundary integrand for the test function gradient term 1926 1927 Note: We are using a first order FEM model for the weak form: 1928 1929 \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 1930 1931 The calling sequence for the callbacks f0 and f1 is given by: 1932 1933 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1934 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1935 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1936 $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]) 1937 1938 + dim - the spatial dimension 1939 . Nf - the number of fields 1940 . uOff - the offset into u[] and u_t[] for each field 1941 . uOff_x - the offset into u_x[] for each field 1942 . u - each field evaluated at the current point 1943 . u_t - the time derivative of each field evaluated at the current point 1944 . u_x - the gradient of each field evaluated at the current point 1945 . aOff - the offset into a[] and a_t[] for each auxiliary field 1946 . aOff_x - the offset into a_x[] for each auxiliary field 1947 . a - each auxiliary field evaluated at the current point 1948 . a_t - the time derivative of each auxiliary field evaluated at the current point 1949 . a_x - the gradient of auxiliary each field evaluated at the current point 1950 . t - current time 1951 . x - coordinates of the current point 1952 . n - unit normal at the current point 1953 . numConstants - number of constant parameters 1954 . constants - constant parameters 1955 - f0 - output values at the current point 1956 1957 Level: intermediate 1958 1959 .seealso: PetscDSSetBdResidual() 1960 @*/ 1961 PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f, 1962 void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1963 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1964 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1965 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), 1966 void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1967 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1968 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1969 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 1970 { 1971 PetscFunctionBegin; 1972 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1973 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 1974 if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];} 1975 if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];} 1976 PetscFunctionReturn(0); 1977 } 1978 1979 /*@C 1980 PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field 1981 1982 Not collective 1983 1984 Input Parameters: 1985 + prob - The PetscDS 1986 . f - The test field number 1987 . f0 - boundary integrand for the test function term 1988 - f1 - boundary integrand for the test function gradient term 1989 1990 Note: We are using a first order FEM model for the weak form: 1991 1992 \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 1993 1994 The calling sequence for the callbacks f0 and f1 is given by: 1995 1996 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1997 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1998 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1999 $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]) 2000 2001 + dim - the spatial dimension 2002 . Nf - the number of fields 2003 . uOff - the offset into u[] and u_t[] for each field 2004 . uOff_x - the offset into u_x[] for each field 2005 . u - each field evaluated at the current point 2006 . u_t - the time derivative of each field evaluated at the current point 2007 . u_x - the gradient of each field evaluated at the current point 2008 . aOff - the offset into a[] and a_t[] for each auxiliary field 2009 . aOff_x - the offset into a_x[] for each auxiliary field 2010 . a - each auxiliary field evaluated at the current point 2011 . a_t - the time derivative of each auxiliary field evaluated at the current point 2012 . a_x - the gradient of auxiliary each field evaluated at the current point 2013 . t - current time 2014 . x - coordinates of the current point 2015 . n - unit normal at the current point 2016 . numConstants - number of constant parameters 2017 . constants - constant parameters 2018 - f0 - output values at the current point 2019 2020 Level: intermediate 2021 2022 .seealso: PetscDSGetBdResidual() 2023 @*/ 2024 PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f, 2025 void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2026 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2027 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2028 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), 2029 void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2030 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2031 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2032 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 2033 { 2034 PetscErrorCode ierr; 2035 2036 PetscFunctionBegin; 2037 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2038 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 2039 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 2040 if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;} 2041 if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;} 2042 PetscFunctionReturn(0); 2043 } 2044 2045 /*@C 2046 PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field 2047 2048 Not collective 2049 2050 Input Parameters: 2051 + prob - The PetscDS 2052 . f - The test field number 2053 - g - The field number 2054 2055 Output Parameters: 2056 + g0 - integrand for the test and basis function term 2057 . g1 - integrand for the test function and basis function gradient term 2058 . g2 - integrand for the test function gradient and basis function term 2059 - g3 - integrand for the test function gradient and basis function gradient term 2060 2061 Note: We are using a first order FEM model for the weak form: 2062 2063 \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 2064 2065 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 2066 2067 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2068 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2069 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2070 $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[]) 2071 2072 + dim - the spatial dimension 2073 . Nf - the number of fields 2074 . uOff - the offset into u[] and u_t[] for each field 2075 . uOff_x - the offset into u_x[] for each field 2076 . u - each field evaluated at the current point 2077 . u_t - the time derivative of each field evaluated at the current point 2078 . u_x - the gradient of each field evaluated at the current point 2079 . aOff - the offset into a[] and a_t[] for each auxiliary field 2080 . aOff_x - the offset into a_x[] for each auxiliary field 2081 . a - each auxiliary field evaluated at the current point 2082 . a_t - the time derivative of each auxiliary field evaluated at the current point 2083 . a_x - the gradient of auxiliary each field evaluated at the current point 2084 . t - current time 2085 . u_tShift - the multiplier a for dF/dU_t 2086 . x - coordinates of the current point 2087 . n - normal at the current point 2088 . numConstants - number of constant parameters 2089 . constants - constant parameters 2090 - g0 - output values at the current point 2091 2092 Level: intermediate 2093 2094 .seealso: PetscDSSetBdJacobian() 2095 @*/ 2096 PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g, 2097 void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2098 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2099 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2100 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 2101 void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2102 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2103 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2104 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 2105 void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2106 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2107 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2108 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 2109 void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2110 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2111 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2112 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 2113 { 2114 PetscFunctionBegin; 2115 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2116 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 2117 if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf); 2118 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];} 2119 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];} 2120 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];} 2121 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];} 2122 PetscFunctionReturn(0); 2123 } 2124 2125 /*@C 2126 PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field 2127 2128 Not collective 2129 2130 Input Parameters: 2131 + prob - The PetscDS 2132 . f - The test field number 2133 . g - The field number 2134 . g0 - integrand for the test and basis function term 2135 . g1 - integrand for the test function and basis function gradient term 2136 . g2 - integrand for the test function gradient and basis function term 2137 - g3 - integrand for the test function gradient and basis function gradient term 2138 2139 Note: We are using a first order FEM model for the weak form: 2140 2141 \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 2142 2143 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 2144 2145 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2146 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2147 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2148 $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[]) 2149 2150 + dim - the spatial dimension 2151 . Nf - the number of fields 2152 . uOff - the offset into u[] and u_t[] for each field 2153 . uOff_x - the offset into u_x[] for each field 2154 . u - each field evaluated at the current point 2155 . u_t - the time derivative of each field evaluated at the current point 2156 . u_x - the gradient of each field evaluated at the current point 2157 . aOff - the offset into a[] and a_t[] for each auxiliary field 2158 . aOff_x - the offset into a_x[] for each auxiliary field 2159 . a - each auxiliary field evaluated at the current point 2160 . a_t - the time derivative of each auxiliary field evaluated at the current point 2161 . a_x - the gradient of auxiliary each field evaluated at the current point 2162 . t - current time 2163 . u_tShift - the multiplier a for dF/dU_t 2164 . x - coordinates of the current point 2165 . n - normal at the current point 2166 . numConstants - number of constant parameters 2167 . constants - constant parameters 2168 - g0 - output values at the current point 2169 2170 Level: intermediate 2171 2172 .seealso: PetscDSGetBdJacobian() 2173 @*/ 2174 PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g, 2175 void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2176 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2177 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2178 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 2179 void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2180 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2181 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2182 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 2183 void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2184 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2185 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2186 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 2187 void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2188 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2189 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2190 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 2191 { 2192 PetscErrorCode ierr; 2193 2194 PetscFunctionBegin; 2195 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2196 if (g0) PetscValidFunction(g0, 4); 2197 if (g1) PetscValidFunction(g1, 5); 2198 if (g2) PetscValidFunction(g2, 6); 2199 if (g3) PetscValidFunction(g3, 7); 2200 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 2201 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 2202 ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr); 2203 prob->gBd[(f*prob->Nf + g)*4+0] = g0; 2204 prob->gBd[(f*prob->Nf + g)*4+1] = g1; 2205 prob->gBd[(f*prob->Nf + g)*4+2] = g2; 2206 prob->gBd[(f*prob->Nf + g)*4+3] = g3; 2207 PetscFunctionReturn(0); 2208 } 2209 2210 /*@C 2211 PetscDSGetExactSolution - Get the pointwise exact solution function for a given test field 2212 2213 Not collective 2214 2215 Input Parameters: 2216 + prob - The PetscDS 2217 - f - The test field number 2218 2219 Output Parameter: 2220 + exactSol - exact solution for the test field 2221 - exactCtx - exact solution context 2222 2223 Note: The calling sequence for the solution functions is given by: 2224 2225 $ sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) 2226 2227 + dim - the spatial dimension 2228 . t - current time 2229 . x - coordinates of the current point 2230 . Nc - the number of field components 2231 . u - the solution field evaluated at the current point 2232 - ctx - a user context 2233 2234 Level: intermediate 2235 2236 .seealso: PetscDSSetExactSolution() 2237 @*/ 2238 PetscErrorCode PetscDSGetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (**sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void **ctx) 2239 { 2240 PetscFunctionBegin; 2241 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2242 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 2243 if (sol) {PetscValidPointer(sol, 3); *sol = prob->exactSol[f];} 2244 if (ctx) {PetscValidPointer(ctx, 4); *ctx = prob->exactCtx[f];} 2245 PetscFunctionReturn(0); 2246 } 2247 2248 /*@C 2249 PetscDSSetExactSolution - Set the pointwise exact solution function for a given test field 2250 2251 Not collective 2252 2253 Input Parameters: 2254 + prob - The PetscDS 2255 . f - The test field number 2256 . sol - solution function for the test fields 2257 - ctx - solution context or NULL 2258 2259 Note: The calling sequence for solution functions is given by: 2260 2261 $ sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) 2262 2263 + dim - the spatial dimension 2264 . t - current time 2265 . x - coordinates of the current point 2266 . Nc - the number of field components 2267 . u - the solution field evaluated at the current point 2268 - ctx - a user context 2269 2270 Level: intermediate 2271 2272 .seealso: PetscDSGetExactSolution() 2273 @*/ 2274 PetscErrorCode PetscDSSetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (*sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void *ctx) 2275 { 2276 PetscErrorCode ierr; 2277 2278 PetscFunctionBegin; 2279 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2280 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 2281 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 2282 if (sol) {PetscValidFunction(sol, 3); prob->exactSol[f] = sol;} 2283 if (ctx) {PetscValidFunction(ctx, 4); prob->exactCtx[f] = ctx;} 2284 PetscFunctionReturn(0); 2285 } 2286 2287 /*@C 2288 PetscDSGetConstants - Returns the array of constants passed to point functions 2289 2290 Not collective 2291 2292 Input Parameter: 2293 . prob - The PetscDS object 2294 2295 Output Parameters: 2296 + numConstants - The number of constants 2297 - constants - The array of constants, NULL if there are none 2298 2299 Level: intermediate 2300 2301 .seealso: PetscDSSetConstants(), PetscDSCreate() 2302 @*/ 2303 PetscErrorCode PetscDSGetConstants(PetscDS prob, PetscInt *numConstants, const PetscScalar *constants[]) 2304 { 2305 PetscFunctionBegin; 2306 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2307 if (numConstants) {PetscValidPointer(numConstants, 2); *numConstants = prob->numConstants;} 2308 if (constants) {PetscValidPointer(constants, 3); *constants = prob->constants;} 2309 PetscFunctionReturn(0); 2310 } 2311 2312 /*@C 2313 PetscDSSetConstants - Set the array of constants passed to point functions 2314 2315 Not collective 2316 2317 Input Parameters: 2318 + prob - The PetscDS object 2319 . numConstants - The number of constants 2320 - constants - The array of constants, NULL if there are none 2321 2322 Level: intermediate 2323 2324 .seealso: PetscDSGetConstants(), PetscDSCreate() 2325 @*/ 2326 PetscErrorCode PetscDSSetConstants(PetscDS prob, PetscInt numConstants, PetscScalar constants[]) 2327 { 2328 PetscErrorCode ierr; 2329 2330 PetscFunctionBegin; 2331 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2332 if (numConstants != prob->numConstants) { 2333 ierr = PetscFree(prob->constants);CHKERRQ(ierr); 2334 prob->numConstants = numConstants; 2335 if (prob->numConstants) { 2336 ierr = PetscMalloc1(prob->numConstants, &prob->constants);CHKERRQ(ierr); 2337 } else { 2338 prob->constants = NULL; 2339 } 2340 } 2341 if (prob->numConstants) { 2342 PetscValidPointer(constants, 3); 2343 ierr = PetscMemcpy(prob->constants, constants, prob->numConstants * sizeof(PetscScalar));CHKERRQ(ierr); 2344 } 2345 PetscFunctionReturn(0); 2346 } 2347 2348 /*@ 2349 PetscDSGetFieldIndex - Returns the index of the given field 2350 2351 Not collective 2352 2353 Input Parameters: 2354 + prob - The PetscDS object 2355 - disc - The discretization object 2356 2357 Output Parameter: 2358 . f - The field number 2359 2360 Level: beginner 2361 2362 .seealso: PetscGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 2363 @*/ 2364 PetscErrorCode PetscDSGetFieldIndex(PetscDS prob, PetscObject disc, PetscInt *f) 2365 { 2366 PetscInt g; 2367 2368 PetscFunctionBegin; 2369 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2370 PetscValidPointer(f, 3); 2371 *f = -1; 2372 for (g = 0; g < prob->Nf; ++g) {if (disc == prob->disc[g]) break;} 2373 if (g == prob->Nf) SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Field not found in PetscDS."); 2374 *f = g; 2375 PetscFunctionReturn(0); 2376 } 2377 2378 /*@ 2379 PetscDSGetFieldSize - Returns the size of the given field in the full space basis 2380 2381 Not collective 2382 2383 Input Parameters: 2384 + prob - The PetscDS object 2385 - f - The field number 2386 2387 Output Parameter: 2388 . size - The size 2389 2390 Level: beginner 2391 2392 .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate() 2393 @*/ 2394 PetscErrorCode PetscDSGetFieldSize(PetscDS prob, PetscInt f, PetscInt *size) 2395 { 2396 PetscErrorCode ierr; 2397 2398 PetscFunctionBegin; 2399 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2400 PetscValidPointer(size, 3); 2401 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 2402 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2403 *size = prob->Nb[f]; 2404 PetscFunctionReturn(0); 2405 } 2406 2407 /*@ 2408 PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis 2409 2410 Not collective 2411 2412 Input Parameters: 2413 + prob - The PetscDS object 2414 - f - The field number 2415 2416 Output Parameter: 2417 . off - The offset 2418 2419 Level: beginner 2420 2421 .seealso: PetscDSGetFieldSize(), PetscDSGetNumFields(), PetscDSCreate() 2422 @*/ 2423 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off) 2424 { 2425 PetscInt size, g; 2426 PetscErrorCode ierr; 2427 2428 PetscFunctionBegin; 2429 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2430 PetscValidPointer(off, 3); 2431 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 2432 *off = 0; 2433 for (g = 0; g < f; ++g) { 2434 ierr = PetscDSGetFieldSize(prob, g, &size);CHKERRQ(ierr); 2435 *off += size; 2436 } 2437 PetscFunctionReturn(0); 2438 } 2439 2440 /*@ 2441 PetscDSGetDimensions - Returns the size of the approximation space for each field on an evaluation point 2442 2443 Not collective 2444 2445 Input Parameter: 2446 . prob - The PetscDS object 2447 2448 Output Parameter: 2449 . dimensions - The number of dimensions 2450 2451 Level: beginner 2452 2453 .seealso: PetscDSGetComponentOffsets(), PetscDSGetNumFields(), PetscDSCreate() 2454 @*/ 2455 PetscErrorCode PetscDSGetDimensions(PetscDS prob, PetscInt *dimensions[]) 2456 { 2457 PetscErrorCode ierr; 2458 2459 PetscFunctionBegin; 2460 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2461 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2462 PetscValidPointer(dimensions, 2); 2463 *dimensions = prob->Nb; 2464 PetscFunctionReturn(0); 2465 } 2466 2467 /*@ 2468 PetscDSGetComponents - Returns the number of components for each field on an evaluation point 2469 2470 Not collective 2471 2472 Input Parameter: 2473 . prob - The PetscDS object 2474 2475 Output Parameter: 2476 . components - The number of components 2477 2478 Level: beginner 2479 2480 .seealso: PetscDSGetComponentOffsets(), PetscDSGetNumFields(), PetscDSCreate() 2481 @*/ 2482 PetscErrorCode PetscDSGetComponents(PetscDS prob, PetscInt *components[]) 2483 { 2484 PetscErrorCode ierr; 2485 2486 PetscFunctionBegin; 2487 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2488 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2489 PetscValidPointer(components, 2); 2490 *components = prob->Nc; 2491 PetscFunctionReturn(0); 2492 } 2493 2494 /*@ 2495 PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point 2496 2497 Not collective 2498 2499 Input Parameters: 2500 + prob - The PetscDS object 2501 - f - The field number 2502 2503 Output Parameter: 2504 . off - The offset 2505 2506 Level: beginner 2507 2508 .seealso: PetscDSGetNumFields(), PetscDSCreate() 2509 @*/ 2510 PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off) 2511 { 2512 PetscFunctionBegin; 2513 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2514 PetscValidPointer(off, 3); 2515 if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf); 2516 *off = prob->off[f]; 2517 PetscFunctionReturn(0); 2518 } 2519 2520 /*@ 2521 PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point 2522 2523 Not collective 2524 2525 Input Parameter: 2526 . prob - The PetscDS object 2527 2528 Output Parameter: 2529 . offsets - The offsets 2530 2531 Level: beginner 2532 2533 .seealso: PetscDSGetNumFields(), PetscDSCreate() 2534 @*/ 2535 PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[]) 2536 { 2537 PetscFunctionBegin; 2538 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2539 PetscValidPointer(offsets, 2); 2540 *offsets = prob->off; 2541 PetscFunctionReturn(0); 2542 } 2543 2544 /*@ 2545 PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point 2546 2547 Not collective 2548 2549 Input Parameter: 2550 . prob - The PetscDS object 2551 2552 Output Parameter: 2553 . offsets - The offsets 2554 2555 Level: beginner 2556 2557 .seealso: PetscDSGetNumFields(), PetscDSCreate() 2558 @*/ 2559 PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[]) 2560 { 2561 PetscFunctionBegin; 2562 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2563 PetscValidPointer(offsets, 2); 2564 *offsets = prob->offDer; 2565 PetscFunctionReturn(0); 2566 } 2567 2568 /*@C 2569 PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization 2570 2571 Not collective 2572 2573 Input Parameter: 2574 . prob - The PetscDS object 2575 2576 Output Parameters: 2577 + basis - The basis function tabulation at quadrature points 2578 - basisDer - The basis function derivative tabulation at quadrature points 2579 2580 Level: intermediate 2581 2582 .seealso: PetscDSCreate() 2583 @*/ 2584 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer) 2585 { 2586 PetscErrorCode ierr; 2587 2588 PetscFunctionBegin; 2589 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2590 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2591 if (basis) {PetscValidPointer(basis, 2); *basis = prob->basis;} 2592 if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;} 2593 PetscFunctionReturn(0); 2594 } 2595 2596 /*@C 2597 PetscDSGetFaceTabulation - Return the basis tabulation at quadrature points on the faces 2598 2599 Not collective 2600 2601 Input Parameter: 2602 . prob - The PetscDS object 2603 2604 Output Parameters: 2605 + basisFace - The basis function tabulation at quadrature points 2606 - basisDerFace - The basis function derivative tabulation at quadrature points 2607 2608 Level: intermediate 2609 2610 .seealso: PetscDSGetTabulation(), PetscDSCreate() 2611 @*/ 2612 PetscErrorCode PetscDSGetFaceTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer) 2613 { 2614 PetscErrorCode ierr; 2615 2616 PetscFunctionBegin; 2617 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2618 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2619 if (basis) {PetscValidPointer(basis, 2); *basis = prob->basisFace;} 2620 if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerFace;} 2621 PetscFunctionReturn(0); 2622 } 2623 2624 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x) 2625 { 2626 PetscErrorCode ierr; 2627 2628 PetscFunctionBegin; 2629 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2630 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2631 if (u) {PetscValidPointer(u, 2); *u = prob->u;} 2632 if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;} 2633 if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;} 2634 PetscFunctionReturn(0); 2635 } 2636 2637 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3) 2638 { 2639 PetscErrorCode ierr; 2640 2641 PetscFunctionBegin; 2642 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2643 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2644 if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;} 2645 if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;} 2646 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;} 2647 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;} 2648 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;} 2649 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;} 2650 PetscFunctionReturn(0); 2651 } 2652 2653 PetscErrorCode PetscDSGetWorkspace(PetscDS prob, PetscReal **x, PetscScalar **basisReal, PetscScalar **basisDerReal, PetscScalar **testReal, PetscScalar **testDerReal) 2654 { 2655 PetscErrorCode ierr; 2656 2657 PetscFunctionBegin; 2658 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2659 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2660 if (x) {PetscValidPointer(x, 2); *x = prob->x;} 2661 if (basisReal) {PetscValidPointer(basisReal, 3); *basisReal = prob->basisReal;} 2662 if (basisDerReal) {PetscValidPointer(basisDerReal, 3); *basisDerReal = prob->basisDerReal;} 2663 if (testReal) {PetscValidPointer(testReal, 3); *testReal = prob->testReal;} 2664 if (testDerReal) {PetscValidPointer(testDerReal, 3); *testDerReal = prob->testDerReal;} 2665 PetscFunctionReturn(0); 2666 } 2667 2668 /*@C 2669 PetscDSAddBoundary - Add a boundary condition to the model 2670 2671 Input Parameters: 2672 + ds - The PetscDS object 2673 . type - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann) 2674 . name - The BC name 2675 . labelname - The label defining constrained points 2676 . field - The field to constrain 2677 . numcomps - The number of constrained field components (0 will constrain all fields) 2678 . comps - An array of constrained component numbers 2679 . bcFunc - A pointwise function giving boundary values 2680 . numids - The number of DMLabel ids for constrained points 2681 . ids - An array of ids for constrained points 2682 - ctx - An optional user context for bcFunc 2683 2684 Options Database Keys: 2685 + -bc_<boundary name> <num> - Overrides the boundary ids 2686 - -bc_<boundary name>_comp <num> - Overrides the boundary components 2687 2688 Level: developer 2689 2690 .seealso: PetscDSGetBoundary() 2691 @*/ 2692 PetscErrorCode PetscDSAddBoundary(PetscDS ds, DMBoundaryConditionType type, const char name[], const char labelname[], PetscInt field, PetscInt numcomps, const PetscInt *comps, void (*bcFunc)(void), PetscInt numids, const PetscInt *ids, void *ctx) 2693 { 2694 DSBoundary b; 2695 PetscErrorCode ierr; 2696 2697 PetscFunctionBegin; 2698 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2699 ierr = PetscNew(&b);CHKERRQ(ierr); 2700 ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr); 2701 ierr = PetscStrallocpy(labelname, (char **) &b->labelname);CHKERRQ(ierr); 2702 ierr = PetscMalloc1(numcomps, &b->comps);CHKERRQ(ierr); 2703 if (numcomps) {ierr = PetscMemcpy(b->comps, comps, numcomps*sizeof(PetscInt));CHKERRQ(ierr);} 2704 ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr); 2705 if (numids) {ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);} 2706 b->type = type; 2707 b->field = field; 2708 b->numcomps = numcomps; 2709 b->func = bcFunc; 2710 b->numids = numids; 2711 b->ctx = ctx; 2712 b->next = ds->boundary; 2713 ds->boundary = b; 2714 PetscFunctionReturn(0); 2715 } 2716 2717 /*@C 2718 PetscDSUpdateBoundary - Change a boundary condition for the model 2719 2720 Input Parameters: 2721 + ds - The PetscDS object 2722 . bd - The boundary condition number 2723 . type - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann) 2724 . name - The BC name 2725 . labelname - The label defining constrained points 2726 . field - The field to constrain 2727 . numcomps - The number of constrained field components 2728 . comps - An array of constrained component numbers 2729 . bcFunc - A pointwise function giving boundary values 2730 . numids - The number of DMLabel ids for constrained points 2731 . ids - An array of ids for constrained points 2732 - ctx - An optional user context for bcFunc 2733 2734 Note: The boundary condition number is the order in which it was registered. The user can get the number of boundary conditions from PetscDSGetNumBoundary(). 2735 2736 Level: developer 2737 2738 .seealso: PetscDSAddBoundary(), PetscDSGetBoundary(), PetscDSGetNumBoundary() 2739 @*/ 2740 PetscErrorCode PetscDSUpdateBoundary(PetscDS ds, PetscInt bd, DMBoundaryConditionType type, const char name[], const char labelname[], PetscInt field, PetscInt numcomps, const PetscInt *comps, void (*bcFunc)(void), PetscInt numids, const PetscInt *ids, void *ctx) 2741 { 2742 DSBoundary b = ds->boundary; 2743 PetscInt n = 0; 2744 PetscErrorCode ierr; 2745 2746 PetscFunctionBegin; 2747 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2748 while (b) { 2749 if (n == bd) break; 2750 b = b->next; 2751 ++n; 2752 } 2753 if (!b) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n); 2754 if (name) { 2755 ierr = PetscFree(b->name);CHKERRQ(ierr); 2756 ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr); 2757 } 2758 if (labelname) { 2759 ierr = PetscFree(b->labelname);CHKERRQ(ierr); 2760 ierr = PetscStrallocpy(labelname, (char **) &b->labelname);CHKERRQ(ierr); 2761 } 2762 if (numcomps >= 0 && numcomps != b->numcomps) { 2763 b->numcomps = numcomps; 2764 ierr = PetscFree(b->comps);CHKERRQ(ierr); 2765 ierr = PetscMalloc1(numcomps, &b->comps);CHKERRQ(ierr); 2766 if (numcomps) {ierr = PetscMemcpy(b->comps, comps, numcomps*sizeof(PetscInt));CHKERRQ(ierr);} 2767 } 2768 if (numids >= 0 && numids != b->numids) { 2769 b->numids = numids; 2770 ierr = PetscFree(b->ids);CHKERRQ(ierr); 2771 ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr); 2772 if (numids) {ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);} 2773 } 2774 b->type = type; 2775 if (field >= 0) {b->field = field;} 2776 if (bcFunc) {b->func = bcFunc;} 2777 if (ctx) {b->ctx = ctx;} 2778 PetscFunctionReturn(0); 2779 } 2780 2781 /*@ 2782 PetscDSGetNumBoundary - Get the number of registered BC 2783 2784 Input Parameters: 2785 . ds - The PetscDS object 2786 2787 Output Parameters: 2788 . numBd - The number of BC 2789 2790 Level: intermediate 2791 2792 .seealso: PetscDSAddBoundary(), PetscDSGetBoundary() 2793 @*/ 2794 PetscErrorCode PetscDSGetNumBoundary(PetscDS ds, PetscInt *numBd) 2795 { 2796 DSBoundary b = ds->boundary; 2797 2798 PetscFunctionBegin; 2799 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2800 PetscValidPointer(numBd, 2); 2801 *numBd = 0; 2802 while (b) {++(*numBd); b = b->next;} 2803 PetscFunctionReturn(0); 2804 } 2805 2806 /*@C 2807 PetscDSGetBoundary - Gets a boundary condition to the model 2808 2809 Input Parameters: 2810 + ds - The PetscDS object 2811 - bd - The BC number 2812 2813 Output Parameters: 2814 + type - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann) 2815 . name - The BC name 2816 . labelname - The label defining constrained points 2817 . field - The field to constrain 2818 . numcomps - The number of constrained field components 2819 . comps - An array of constrained component numbers 2820 . bcFunc - A pointwise function giving boundary values 2821 . numids - The number of DMLabel ids for constrained points 2822 . ids - An array of ids for constrained points 2823 - ctx - An optional user context for bcFunc 2824 2825 Options Database Keys: 2826 + -bc_<boundary name> <num> - Overrides the boundary ids 2827 - -bc_<boundary name>_comp <num> - Overrides the boundary components 2828 2829 Level: developer 2830 2831 .seealso: PetscDSAddBoundary() 2832 @*/ 2833 PetscErrorCode PetscDSGetBoundary(PetscDS ds, PetscInt bd, DMBoundaryConditionType *type, const char **name, const char **labelname, PetscInt *field, PetscInt *numcomps, const PetscInt **comps, void (**func)(void), PetscInt *numids, const PetscInt **ids, void **ctx) 2834 { 2835 DSBoundary b = ds->boundary; 2836 PetscInt n = 0; 2837 2838 PetscFunctionBegin; 2839 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2840 while (b) { 2841 if (n == bd) break; 2842 b = b->next; 2843 ++n; 2844 } 2845 if (!b) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n); 2846 if (type) { 2847 PetscValidPointer(type, 3); 2848 *type = b->type; 2849 } 2850 if (name) { 2851 PetscValidPointer(name, 4); 2852 *name = b->name; 2853 } 2854 if (labelname) { 2855 PetscValidPointer(labelname, 5); 2856 *labelname = b->labelname; 2857 } 2858 if (field) { 2859 PetscValidPointer(field, 6); 2860 *field = b->field; 2861 } 2862 if (numcomps) { 2863 PetscValidPointer(numcomps, 7); 2864 *numcomps = b->numcomps; 2865 } 2866 if (comps) { 2867 PetscValidPointer(comps, 8); 2868 *comps = b->comps; 2869 } 2870 if (func) { 2871 PetscValidPointer(func, 9); 2872 *func = b->func; 2873 } 2874 if (numids) { 2875 PetscValidPointer(numids, 10); 2876 *numids = b->numids; 2877 } 2878 if (ids) { 2879 PetscValidPointer(ids, 11); 2880 *ids = b->ids; 2881 } 2882 if (ctx) { 2883 PetscValidPointer(ctx, 12); 2884 *ctx = b->ctx; 2885 } 2886 PetscFunctionReturn(0); 2887 } 2888 2889 /*@ 2890 PetscDSCopyBoundary - Copy all boundary condition objects to the new problem 2891 2892 Not collective 2893 2894 Input Parameter: 2895 . prob - The PetscDS object 2896 2897 Output Parameter: 2898 . newprob - The PetscDS copy 2899 2900 Level: intermediate 2901 2902 .seealso: PetscDSCopyEquations(), PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate() 2903 @*/ 2904 PetscErrorCode PetscDSCopyBoundary(PetscDS probA, PetscDS probB) 2905 { 2906 DSBoundary b, next, *lastnext; 2907 PetscErrorCode ierr; 2908 2909 PetscFunctionBegin; 2910 PetscValidHeaderSpecific(probA, PETSCDS_CLASSID, 1); 2911 PetscValidHeaderSpecific(probB, PETSCDS_CLASSID, 2); 2912 if (probA == probB) PetscFunctionReturn(0); 2913 next = probB->boundary; 2914 while (next) { 2915 DSBoundary b = next; 2916 2917 next = b->next; 2918 ierr = PetscFree(b->comps);CHKERRQ(ierr); 2919 ierr = PetscFree(b->ids);CHKERRQ(ierr); 2920 ierr = PetscFree(b->name);CHKERRQ(ierr); 2921 ierr = PetscFree(b->labelname);CHKERRQ(ierr); 2922 ierr = PetscFree(b);CHKERRQ(ierr); 2923 } 2924 lastnext = &(probB->boundary); 2925 for (b = probA->boundary; b; b = b->next) { 2926 DSBoundary bNew; 2927 2928 ierr = PetscNew(&bNew);CHKERRQ(ierr); 2929 bNew->numcomps = b->numcomps; 2930 ierr = PetscMalloc1(bNew->numcomps, &bNew->comps);CHKERRQ(ierr); 2931 ierr = PetscMemcpy(bNew->comps, b->comps, bNew->numcomps*sizeof(PetscInt));CHKERRQ(ierr); 2932 bNew->numids = b->numids; 2933 ierr = PetscMalloc1(bNew->numids, &bNew->ids);CHKERRQ(ierr); 2934 ierr = PetscMemcpy(bNew->ids, b->ids, bNew->numids*sizeof(PetscInt));CHKERRQ(ierr); 2935 ierr = PetscStrallocpy(b->labelname,(char **) &(bNew->labelname));CHKERRQ(ierr); 2936 ierr = PetscStrallocpy(b->name,(char **) &(bNew->name));CHKERRQ(ierr); 2937 bNew->ctx = b->ctx; 2938 bNew->type = b->type; 2939 bNew->field = b->field; 2940 bNew->func = b->func; 2941 2942 *lastnext = bNew; 2943 lastnext = &(bNew->next); 2944 } 2945 PetscFunctionReturn(0); 2946 } 2947 2948 /*@C 2949 PetscDSSelectEquations - Copy pointwise function pointers to the new problem with different field layout 2950 2951 Not collective 2952 2953 Input Parameter: 2954 + prob - The PetscDS object 2955 . numFields - Number of new fields 2956 - fields - Old field number for each new field 2957 2958 Output Parameter: 2959 . newprob - The PetscDS copy 2960 2961 Level: intermediate 2962 2963 .seealso: PetscDSCopyBoundary(), PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate() 2964 @*/ 2965 PetscErrorCode PetscDSSelectEquations(PetscDS prob, PetscInt numFields, const PetscInt fields[], PetscDS newprob) 2966 { 2967 PetscInt Nf, Nfn, fn, gn; 2968 PetscErrorCode ierr; 2969 2970 PetscFunctionBegin; 2971 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2972 if (fields) PetscValidPointer(fields, 3); 2973 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 4); 2974 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 2975 ierr = PetscDSGetNumFields(newprob, &Nfn);CHKERRQ(ierr); 2976 if (numFields > Nfn) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_SIZ, "Number of fields %D to transfer must not be greater then the total number of fields %D", numFields, Nfn); 2977 for (fn = 0; fn < numFields; ++fn) { 2978 const PetscInt f = fields ? fields[fn] : fn; 2979 PetscPointFunc obj; 2980 PetscPointFunc f0, f1; 2981 PetscBdPointFunc f0Bd, f1Bd; 2982 PetscRiemannFunc r; 2983 2984 if (f >= Nf) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_SIZ, "Field %D must be in [0, %D)", f, Nf); 2985 ierr = PetscDSGetObjective(prob, f, &obj);CHKERRQ(ierr); 2986 ierr = PetscDSGetResidual(prob, f, &f0, &f1);CHKERRQ(ierr); 2987 ierr = PetscDSGetBdResidual(prob, f, &f0Bd, &f1Bd);CHKERRQ(ierr); 2988 ierr = PetscDSGetRiemannSolver(prob, f, &r);CHKERRQ(ierr); 2989 ierr = PetscDSSetObjective(newprob, fn, obj);CHKERRQ(ierr); 2990 ierr = PetscDSSetResidual(newprob, fn, f0, f1);CHKERRQ(ierr); 2991 ierr = PetscDSSetBdResidual(newprob, fn, f0Bd, f1Bd);CHKERRQ(ierr); 2992 ierr = PetscDSSetRiemannSolver(newprob, fn, r);CHKERRQ(ierr); 2993 for (gn = 0; gn < numFields; ++gn) { 2994 const PetscInt g = fields ? fields[gn] : gn; 2995 PetscPointJac g0, g1, g2, g3; 2996 PetscPointJac g0p, g1p, g2p, g3p; 2997 PetscBdPointJac g0Bd, g1Bd, g2Bd, g3Bd; 2998 2999 if (g >= Nf) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_SIZ, "Field %D must be in [0, %D)", g, Nf); 3000 ierr = PetscDSGetJacobian(prob, f, g, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 3001 ierr = PetscDSGetJacobianPreconditioner(prob, f, g, &g0p, &g1p, &g2p, &g3p);CHKERRQ(ierr); 3002 ierr = PetscDSGetBdJacobian(prob, f, g, &g0Bd, &g1Bd, &g2Bd, &g3Bd);CHKERRQ(ierr); 3003 ierr = PetscDSSetJacobian(newprob, fn, gn, g0, g1, g2, g3);CHKERRQ(ierr); 3004 ierr = PetscDSSetJacobianPreconditioner(prob, fn, gn, g0p, g1p, g2p, g3p);CHKERRQ(ierr); 3005 ierr = PetscDSSetBdJacobian(newprob, fn, gn, g0Bd, g1Bd, g2Bd, g3Bd);CHKERRQ(ierr); 3006 } 3007 } 3008 PetscFunctionReturn(0); 3009 } 3010 3011 /*@ 3012 PetscDSCopyEquations - Copy all pointwise function pointers to the new problem 3013 3014 Not collective 3015 3016 Input Parameter: 3017 . prob - The PetscDS object 3018 3019 Output Parameter: 3020 . newprob - The PetscDS copy 3021 3022 Level: intermediate 3023 3024 .seealso: PetscDSCopyBoundary(), PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate() 3025 @*/ 3026 PetscErrorCode PetscDSCopyEquations(PetscDS prob, PetscDS newprob) 3027 { 3028 PetscInt Nf, Ng; 3029 PetscErrorCode ierr; 3030 3031 PetscFunctionBegin; 3032 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3033 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2); 3034 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 3035 ierr = PetscDSGetNumFields(newprob, &Ng);CHKERRQ(ierr); 3036 if (Nf != Ng) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_SIZ, "Number of fields must match %D != %D", Nf, Ng); 3037 ierr = PetscDSSelectEquations(prob, Nf, NULL, newprob);CHKERRQ(ierr); 3038 PetscFunctionReturn(0); 3039 } 3040 /*@ 3041 PetscDSCopyConstants - Copy all constants to the new problem 3042 3043 Not collective 3044 3045 Input Parameter: 3046 . prob - The PetscDS object 3047 3048 Output Parameter: 3049 . newprob - The PetscDS copy 3050 3051 Level: intermediate 3052 3053 .seealso: PetscDSCopyBoundary(), PetscDSCopyEquations(), PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate() 3054 @*/ 3055 PetscErrorCode PetscDSCopyConstants(PetscDS prob, PetscDS newprob) 3056 { 3057 PetscInt Nc; 3058 const PetscScalar *constants; 3059 PetscErrorCode ierr; 3060 3061 PetscFunctionBegin; 3062 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3063 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2); 3064 ierr = PetscDSGetConstants(prob, &Nc, &constants);CHKERRQ(ierr); 3065 ierr = PetscDSSetConstants(newprob, Nc, (PetscScalar *) constants);CHKERRQ(ierr); 3066 PetscFunctionReturn(0); 3067 } 3068 3069 PetscErrorCode PetscDSGetHeightSubspace(PetscDS prob, PetscInt height, PetscDS *subprob) 3070 { 3071 PetscInt dim, Nf, f; 3072 PetscErrorCode ierr; 3073 3074 PetscFunctionBegin; 3075 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3076 PetscValidPointer(subprob, 3); 3077 if (height == 0) {*subprob = prob; PetscFunctionReturn(0);} 3078 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 3079 ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr); 3080 if (height > dim) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_OUTOFRANGE, "DS can only handle height in [0, %D], not %D", dim, height); 3081 if (!prob->subprobs) {ierr = PetscCalloc1(dim, &prob->subprobs);CHKERRQ(ierr);} 3082 if (!prob->subprobs[height-1]) { 3083 PetscInt cdim; 3084 3085 ierr = PetscDSCreate(PetscObjectComm((PetscObject) prob), &prob->subprobs[height-1]);CHKERRQ(ierr); 3086 ierr = PetscDSGetCoordinateDimension(prob, &cdim);CHKERRQ(ierr); 3087 ierr = PetscDSSetCoordinateDimension(prob->subprobs[height-1], cdim);CHKERRQ(ierr); 3088 for (f = 0; f < Nf; ++f) { 3089 PetscFE subfe; 3090 PetscObject obj; 3091 PetscClassId id; 3092 3093 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 3094 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 3095 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetHeightSubspace((PetscFE) obj, height, &subfe);CHKERRQ(ierr);} 3096 else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported discretization type for field %d", f); 3097 ierr = PetscDSSetDiscretization(prob->subprobs[height-1], f, (PetscObject) subfe);CHKERRQ(ierr); 3098 } 3099 } 3100 *subprob = prob->subprobs[height-1]; 3101 PetscFunctionReturn(0); 3102 } 3103 3104 static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob) 3105 { 3106 PetscErrorCode ierr; 3107 3108 PetscFunctionBegin; 3109 ierr = PetscFree(prob->data);CHKERRQ(ierr); 3110 PetscFunctionReturn(0); 3111 } 3112 3113 static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob) 3114 { 3115 PetscFunctionBegin; 3116 prob->ops->setfromoptions = NULL; 3117 prob->ops->setup = NULL; 3118 prob->ops->view = NULL; 3119 prob->ops->destroy = PetscDSDestroy_Basic; 3120 PetscFunctionReturn(0); 3121 } 3122 3123 /*MC 3124 PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions 3125 3126 Level: intermediate 3127 3128 .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType() 3129 M*/ 3130 3131 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob) 3132 { 3133 PetscDS_Basic *b; 3134 PetscErrorCode ierr; 3135 3136 PetscFunctionBegin; 3137 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3138 ierr = PetscNewLog(prob, &b);CHKERRQ(ierr); 3139 prob->data = b; 3140 3141 ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr); 3142 PetscFunctionReturn(0); 3143 } 3144