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