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