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