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