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