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