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