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