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