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