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