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