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