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 if (!PetscDSRegisterAllCalled) {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 if (!PetscDSRegisterAllCalled) {ierr = PetscDSRegisterAll();CHKERRQ(ierr);} 121 *name = ((PetscObject) prob)->type_name; 122 PetscFunctionReturn(0); 123 } 124 125 #undef __FUNCT__ 126 #define __FUNCT__ "PetscDSView" 127 /*@C 128 PetscDSView - Views a PetscDS 129 130 Collective on PetscDS 131 132 Input Parameter: 133 + prob - the PetscDS object to view 134 - v - the viewer 135 136 Level: developer 137 138 .seealso PetscDSDestroy() 139 @*/ 140 PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v) 141 { 142 PetscErrorCode ierr; 143 144 PetscFunctionBegin; 145 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 146 if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);} 147 if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);} 148 PetscFunctionReturn(0); 149 } 150 151 #undef __FUNCT__ 152 #define __FUNCT__ "PetscDSViewFromOptions" 153 /* 154 PetscDSViewFromOptions - Processes command line options to determine if/how a PetscDS is to be viewed. 155 156 Collective on PetscDS 157 158 Input Parameters: 159 + prob - the PetscDS 160 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd' 161 - optionname - option to activate viewing 162 163 Level: intermediate 164 165 .keywords: PetscDS, view, options, database 166 .seealso: VecViewFromOptions(), MatViewFromOptions() 167 */ 168 PetscErrorCode PetscDSViewFromOptions(PetscDS prob, const char prefix[], const char optionname[]) 169 { 170 PetscViewer viewer; 171 PetscViewerFormat format; 172 PetscBool flg; 173 PetscErrorCode ierr; 174 175 PetscFunctionBegin; 176 if (prefix) {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);} 177 else {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), ((PetscObject) prob)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);} 178 if (flg) { 179 ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 180 ierr = PetscDSView(prob, viewer);CHKERRQ(ierr); 181 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 182 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 183 } 184 PetscFunctionReturn(0); 185 } 186 187 #undef __FUNCT__ 188 #define __FUNCT__ "PetscDSSetFromOptions" 189 /*@ 190 PetscDSSetFromOptions - sets parameters in a PetscDS from the options database 191 192 Collective on PetscDS 193 194 Input Parameter: 195 . prob - the PetscDS object to set options for 196 197 Options Database: 198 199 Level: developer 200 201 .seealso PetscDSView() 202 @*/ 203 PetscErrorCode PetscDSSetFromOptions(PetscDS prob) 204 { 205 const char *defaultType; 206 char name[256]; 207 PetscBool flg; 208 PetscErrorCode ierr; 209 210 PetscFunctionBegin; 211 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 212 if (!((PetscObject) prob)->type_name) { 213 defaultType = PETSCDSBASIC; 214 } else { 215 defaultType = ((PetscObject) prob)->type_name; 216 } 217 if (!PetscDSRegisterAllCalled) {ierr = PetscDSRegisterAll();CHKERRQ(ierr);} 218 219 ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr); 220 ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr); 221 if (flg) { 222 ierr = PetscDSSetType(prob, name);CHKERRQ(ierr); 223 } else if (!((PetscObject) prob)->type_name) { 224 ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr); 225 } 226 if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);} 227 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 228 ierr = PetscObjectProcessOptionsHandlers((PetscObject) prob);CHKERRQ(ierr); 229 ierr = PetscOptionsEnd();CHKERRQ(ierr); 230 ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr); 231 PetscFunctionReturn(0); 232 } 233 234 #undef __FUNCT__ 235 #define __FUNCT__ "PetscDSSetUp" 236 /*@C 237 PetscDSSetUp - Construct data structures for the PetscDS 238 239 Collective on PetscDS 240 241 Input Parameter: 242 . prob - the PetscDS object to setup 243 244 Level: developer 245 246 .seealso PetscDSView(), PetscDSDestroy() 247 @*/ 248 PetscErrorCode PetscDSSetUp(PetscDS prob) 249 { 250 const PetscInt Nf = prob->Nf; 251 PetscInt dim, work, NcMax = 0, NqMax = 0, f; 252 PetscErrorCode ierr; 253 254 PetscFunctionBegin; 255 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 256 if (prob->setup) PetscFunctionReturn(0); 257 /* Calculate sizes */ 258 ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr); 259 prob->totDim = prob->totDimBd = prob->totComp = 0; 260 ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisBd,Nf,&prob->basisDerBd);CHKERRQ(ierr); 261 for (f = 0; f < Nf; ++f) { 262 PetscFE fe = (PetscFE) prob->disc[f]; 263 PetscFE feBd = (PetscFE) prob->discBd[f]; 264 PetscQuadrature q; 265 PetscInt Nq, Nb, Nc; 266 267 /* TODO Dispatch on discretization type*/ 268 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 269 ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 270 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 271 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 272 ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr); 273 NqMax = PetscMax(NqMax, Nq); 274 NcMax = PetscMax(NcMax, Nc); 275 prob->totDim += Nb*Nc; 276 prob->totComp += Nc; 277 if (feBd) { 278 ierr = PetscFEGetDimension(feBd, &Nb);CHKERRQ(ierr); 279 ierr = PetscFEGetNumComponents(feBd, &Nc);CHKERRQ(ierr); 280 ierr = PetscFEGetDefaultTabulation(feBd, &prob->basisBd[f], &prob->basisDerBd[f], NULL);CHKERRQ(ierr); 281 prob->totDimBd += Nb*Nc; 282 } 283 } 284 work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim)); 285 /* Allocate works space */ 286 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); 287 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); 288 if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);} 289 prob->setup = PETSC_TRUE; 290 PetscFunctionReturn(0); 291 } 292 293 #undef __FUNCT__ 294 #define __FUNCT__ "PetscDSDestroyStructs_Static" 295 static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob) 296 { 297 PetscErrorCode ierr; 298 299 PetscFunctionBegin; 300 ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisBd,prob->basisDerBd);CHKERRQ(ierr); 301 ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr); 302 ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr); 303 PetscFunctionReturn(0); 304 } 305 306 #undef __FUNCT__ 307 #define __FUNCT__ "PetscDSEnlarge_Static" 308 static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew) 309 { 310 PetscObject *tmpd, *tmpdbd; 311 PointFunc *tmpobj, *tmpf, *tmpg; 312 BdPointFunc *tmpfbd, *tmpgbd; 313 RiemannFunc *tmpr; 314 void **tmpctx; 315 PetscInt Nf = prob->Nf, f; 316 PetscErrorCode ierr; 317 318 PetscFunctionBegin; 319 if (Nf >= NfNew) PetscFunctionReturn(0); 320 prob->setup = PETSC_FALSE; 321 ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr); 322 ierr = PetscMalloc1(NfNew, &tmpd);CHKERRQ(ierr); 323 for (f = 0; f < Nf; ++f) tmpd[f] = prob->disc[f]; 324 for (f = Nf; f < NfNew; ++f) {ierr = PetscContainerCreate(PetscObjectComm((PetscObject) prob), (PetscContainer *) &tmpd[f]);CHKERRQ(ierr);} 325 ierr = PetscFree(prob->disc);CHKERRQ(ierr); 326 prob->Nf = NfNew; 327 prob->disc = tmpd; 328 ierr = PetscCalloc5(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr); 329 for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f]; 330 for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f]; 331 for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f]; 332 for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f]; 333 for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f]; 334 for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL; 335 for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL; 336 for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL; 337 for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL; 338 for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL; 339 ierr = PetscFree5(prob->obj, prob->f, prob->g, prob->r, prob->ctx);CHKERRQ(ierr); 340 prob->obj = tmpobj; 341 prob->f = tmpf; 342 prob->g = tmpg; 343 prob->r = tmpr; 344 prob->ctx = tmpctx; 345 ierr = PetscMalloc1(NfNew, &tmpdbd);CHKERRQ(ierr); 346 for (f = 0; f < Nf; ++f) tmpdbd[f] = prob->discBd[f]; 347 for (f = Nf; f < NfNew; ++f) tmpdbd[f] = NULL; 348 ierr = PetscFree(prob->discBd);CHKERRQ(ierr); 349 prob->discBd = tmpdbd; 350 ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr); 351 for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f]; 352 for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f]; 353 for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL; 354 for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL; 355 ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr); 356 prob->fBd = tmpfbd; 357 prob->gBd = tmpgbd; 358 PetscFunctionReturn(0); 359 } 360 361 #undef __FUNCT__ 362 #define __FUNCT__ "PetscDSDestroy" 363 /*@ 364 PetscDSDestroy - Destroys a PetscDS object 365 366 Collective on PetscDS 367 368 Input Parameter: 369 . prob - the PetscDS object to destroy 370 371 Level: developer 372 373 .seealso PetscDSView() 374 @*/ 375 PetscErrorCode PetscDSDestroy(PetscDS *prob) 376 { 377 PetscInt f; 378 PetscErrorCode ierr; 379 380 PetscFunctionBegin; 381 if (!*prob) PetscFunctionReturn(0); 382 PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1); 383 384 if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);} 385 ((PetscObject) (*prob))->refct = 0; 386 ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr); 387 for (f = 0; f < (*prob)->Nf; ++f) { 388 ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr); 389 ierr = PetscObjectDereference((*prob)->discBd[f]);CHKERRQ(ierr); 390 } 391 ierr = PetscFree((*prob)->disc);CHKERRQ(ierr); 392 ierr = PetscFree((*prob)->discBd);CHKERRQ(ierr); 393 ierr = PetscFree5((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr); 394 ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr); 395 if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);} 396 ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr); 397 PetscFunctionReturn(0); 398 } 399 400 #undef __FUNCT__ 401 #define __FUNCT__ "PetscDSCreate" 402 /*@ 403 PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType(). 404 405 Collective on MPI_Comm 406 407 Input Parameter: 408 . comm - The communicator for the PetscDS object 409 410 Output Parameter: 411 . prob - The PetscDS object 412 413 Level: beginner 414 415 .seealso: PetscDSSetType(), PETSCDSBASIC 416 @*/ 417 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob) 418 { 419 PetscDS p; 420 PetscErrorCode ierr; 421 422 PetscFunctionBegin; 423 PetscValidPointer(prob, 2); 424 *prob = NULL; 425 ierr = PetscDSInitializePackage();CHKERRQ(ierr); 426 427 ierr = PetscHeaderCreate(p, _p_PetscDS, struct _PetscDSOps, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr); 428 ierr = PetscMemzero(p->ops, sizeof(struct _PetscDSOps));CHKERRQ(ierr); 429 430 p->Nf = 0; 431 p->setup = PETSC_FALSE; 432 433 *prob = p; 434 PetscFunctionReturn(0); 435 } 436 437 #undef __FUNCT__ 438 #define __FUNCT__ "PetscDSGetNumFields" 439 PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf) 440 { 441 PetscFunctionBegin; 442 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 443 PetscValidPointer(Nf, 2); 444 *Nf = prob->Nf; 445 PetscFunctionReturn(0); 446 } 447 448 #undef __FUNCT__ 449 #define __FUNCT__ "PetscDSGetSpatialDimension" 450 PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim) 451 { 452 PetscErrorCode ierr; 453 454 PetscFunctionBegin; 455 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 456 PetscValidPointer(dim, 2); 457 *dim = 0; 458 if (prob->Nf) {ierr = PetscFEGetSpatialDimension((PetscFE) prob->disc[0], dim);CHKERRQ(ierr);} 459 PetscFunctionReturn(0); 460 } 461 462 #undef __FUNCT__ 463 #define __FUNCT__ "PetscDSGetTotalDimension" 464 PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim) 465 { 466 PetscErrorCode ierr; 467 468 PetscFunctionBegin; 469 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 470 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 471 PetscValidPointer(dim, 2); 472 *dim = prob->totDim; 473 PetscFunctionReturn(0); 474 } 475 476 #undef __FUNCT__ 477 #define __FUNCT__ "PetscDSGetTotalBdDimension" 478 PetscErrorCode PetscDSGetTotalBdDimension(PetscDS prob, PetscInt *dim) 479 { 480 PetscErrorCode ierr; 481 482 PetscFunctionBegin; 483 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 484 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 485 PetscValidPointer(dim, 2); 486 *dim = prob->totDimBd; 487 PetscFunctionReturn(0); 488 } 489 490 #undef __FUNCT__ 491 #define __FUNCT__ "PetscDSGetTotalComponents" 492 PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc) 493 { 494 PetscErrorCode ierr; 495 496 PetscFunctionBegin; 497 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 498 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 499 PetscValidPointer(Nc, 2); 500 *Nc = prob->totComp; 501 PetscFunctionReturn(0); 502 } 503 504 #undef __FUNCT__ 505 #define __FUNCT__ "PetscDSGetDiscretization" 506 PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc) 507 { 508 PetscFunctionBegin; 509 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 510 PetscValidPointer(disc, 3); 511 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); 512 *disc = prob->disc[f]; 513 PetscFunctionReturn(0); 514 } 515 516 #undef __FUNCT__ 517 #define __FUNCT__ "PetscDSGetBdDiscretization" 518 PetscErrorCode PetscDSGetBdDiscretization(PetscDS prob, PetscInt f, PetscObject *disc) 519 { 520 PetscFunctionBegin; 521 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 522 PetscValidPointer(disc, 3); 523 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); 524 *disc = prob->discBd[f]; 525 PetscFunctionReturn(0); 526 } 527 528 #undef __FUNCT__ 529 #define __FUNCT__ "PetscDSSetDiscretization" 530 PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc) 531 { 532 PetscErrorCode ierr; 533 534 PetscFunctionBegin; 535 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 536 PetscValidPointer(disc, 3); 537 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 538 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 539 if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);} 540 prob->disc[f] = disc; 541 ierr = PetscObjectReference(disc);CHKERRQ(ierr); 542 PetscFunctionReturn(0); 543 } 544 545 #undef __FUNCT__ 546 #define __FUNCT__ "PetscDSSetBdDiscretization" 547 PetscErrorCode PetscDSSetBdDiscretization(PetscDS prob, PetscInt f, PetscObject disc) 548 { 549 PetscErrorCode ierr; 550 551 PetscFunctionBegin; 552 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 553 if (disc) PetscValidPointer(disc, 3); 554 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 555 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 556 if (prob->discBd[f]) {ierr = PetscObjectDereference(prob->discBd[f]);CHKERRQ(ierr);} 557 prob->discBd[f] = disc; 558 ierr = PetscObjectReference(disc);CHKERRQ(ierr); 559 PetscFunctionReturn(0); 560 } 561 562 #undef __FUNCT__ 563 #define __FUNCT__ "PetscDSAddDiscretization" 564 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc) 565 { 566 PetscErrorCode ierr; 567 568 PetscFunctionBegin; 569 ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr); 570 PetscFunctionReturn(0); 571 } 572 573 #undef __FUNCT__ 574 #define __FUNCT__ "PetscDSAddBdDiscretization" 575 PetscErrorCode PetscDSAddBdDiscretization(PetscDS prob, PetscObject disc) 576 { 577 PetscErrorCode ierr; 578 579 PetscFunctionBegin; 580 ierr = PetscDSSetBdDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr); 581 PetscFunctionReturn(0); 582 } 583 584 #undef __FUNCT__ 585 #define __FUNCT__ "PetscDSGetObjective" 586 PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f, 587 void (**obj)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar obj[])) 588 { 589 PetscFunctionBegin; 590 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 591 PetscValidPointer(obj, 2); 592 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); 593 *obj = prob->obj[f]; 594 PetscFunctionReturn(0); 595 } 596 597 #undef __FUNCT__ 598 #define __FUNCT__ "PetscDSSetObjective" 599 PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f, 600 void (*obj)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar obj[])) 601 { 602 PetscErrorCode ierr; 603 604 PetscFunctionBegin; 605 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 606 PetscValidFunction(obj, 2); 607 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 608 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 609 prob->obj[f] = obj; 610 PetscFunctionReturn(0); 611 } 612 613 #undef __FUNCT__ 614 #define __FUNCT__ "PetscDSGetResidual" 615 PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f, 616 void (**f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f0[]), 617 void (**f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f1[])) 618 { 619 PetscFunctionBegin; 620 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 621 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); 622 if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];} 623 if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];} 624 PetscFunctionReturn(0); 625 } 626 627 #undef __FUNCT__ 628 #define __FUNCT__ "PetscDSSetResidual" 629 PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f, 630 void (*f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f0[]), 631 void (*f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f1[])) 632 { 633 PetscErrorCode ierr; 634 635 PetscFunctionBegin; 636 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 637 PetscValidFunction(f0, 3); 638 PetscValidFunction(f1, 4); 639 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 640 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 641 prob->f[f*2+0] = f0; 642 prob->f[f*2+1] = f1; 643 PetscFunctionReturn(0); 644 } 645 646 #undef __FUNCT__ 647 #define __FUNCT__ "PetscDSGetJacobian" 648 PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g, 649 void (**g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g0[]), 650 void (**g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g1[]), 651 void (**g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g2[]), 652 void (**g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g3[])) 653 { 654 PetscFunctionBegin; 655 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 656 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); 657 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); 658 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];} 659 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];} 660 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];} 661 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];} 662 PetscFunctionReturn(0); 663 } 664 665 #undef __FUNCT__ 666 #define __FUNCT__ "PetscDSSetJacobian" 667 PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g, 668 void (*g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g0[]), 669 void (*g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g1[]), 670 void (*g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g2[]), 671 void (*g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g3[])) 672 { 673 PetscErrorCode ierr; 674 675 PetscFunctionBegin; 676 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 677 if (g0) PetscValidFunction(g0, 4); 678 if (g1) PetscValidFunction(g1, 5); 679 if (g2) PetscValidFunction(g2, 6); 680 if (g3) PetscValidFunction(g3, 7); 681 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 682 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 683 ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr); 684 prob->g[(f*prob->Nf + g)*4+0] = g0; 685 prob->g[(f*prob->Nf + g)*4+1] = g1; 686 prob->g[(f*prob->Nf + g)*4+2] = g2; 687 prob->g[(f*prob->Nf + g)*4+3] = g3; 688 PetscFunctionReturn(0); 689 } 690 691 #undef __FUNCT__ 692 #define __FUNCT__ "PetscDSGetRiemannSolver" 693 /*@C 694 PetscDSGetRiemannSolver - Returns the Riemann solver for the given field 695 696 Not collective 697 698 Input Arguments: 699 + prob - The PetscDS object 700 - f - The field number 701 702 Output Argument: 703 . r - Riemann solver 704 705 Calling sequence for r: 706 707 $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 708 709 + x - The coordinates at a point on the interface 710 . n - The normal vector to the interface 711 . uL - The state vector to the left of the interface 712 . uR - The state vector to the right of the interface 713 . flux - output array of flux through the interface 714 - ctx - optional user context 715 716 Level: intermediate 717 718 .seealso: PetscDSSetRiemannSolver() 719 @*/ 720 PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f, 721 void (**r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)) 722 { 723 PetscFunctionBegin; 724 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 725 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); 726 PetscValidPointer(r, 3); 727 *r = prob->r[f]; 728 PetscFunctionReturn(0); 729 } 730 731 #undef __FUNCT__ 732 #define __FUNCT__ "PetscDSSetRiemannSolver" 733 /*@C 734 PetscDSSetRiemannSolver - Sets the Riemann solver for the given field 735 736 Not collective 737 738 Input Arguments: 739 + prob - The PetscDS object 740 . f - The field number 741 - r - Riemann solver 742 743 Calling sequence for r: 744 745 $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 746 747 + x - The coordinates at a point on the interface 748 . n - The normal vector to the interface 749 . uL - The state vector to the left of the interface 750 . uR - The state vector to the right of the interface 751 . flux - output array of flux through the interface 752 - ctx - optional user context 753 754 Level: intermediate 755 756 .seealso: PetscDSGetRiemannSolver() 757 @*/ 758 PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f, 759 void (*r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)) 760 { 761 PetscErrorCode ierr; 762 763 PetscFunctionBegin; 764 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 765 PetscValidFunction(r, 3); 766 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 767 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 768 prob->r[f] = r; 769 PetscFunctionReturn(0); 770 } 771 772 #undef __FUNCT__ 773 #define __FUNCT__ "PetscDSGetContext" 774 PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx) 775 { 776 PetscFunctionBegin; 777 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 778 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); 779 PetscValidPointer(ctx, 3); 780 *ctx = prob->ctx[f]; 781 PetscFunctionReturn(0); 782 } 783 784 #undef __FUNCT__ 785 #define __FUNCT__ "PetscDSSetContext" 786 PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx) 787 { 788 PetscErrorCode ierr; 789 790 PetscFunctionBegin; 791 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 792 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 793 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 794 prob->ctx[f] = ctx; 795 PetscFunctionReturn(0); 796 } 797 798 #undef __FUNCT__ 799 #define __FUNCT__ "PetscDSGetBdResidual" 800 PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f, 801 void (**f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f0[]), 802 void (**f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f1[])) 803 { 804 PetscFunctionBegin; 805 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 806 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); 807 if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];} 808 if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];} 809 PetscFunctionReturn(0); 810 } 811 812 #undef __FUNCT__ 813 #define __FUNCT__ "PetscDSSetBdResidual" 814 PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f, 815 void (*f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f0[]), 816 void (*f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f1[])) 817 { 818 PetscErrorCode ierr; 819 820 PetscFunctionBegin; 821 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 822 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 823 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 824 if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;} 825 if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;} 826 PetscFunctionReturn(0); 827 } 828 829 #undef __FUNCT__ 830 #define __FUNCT__ "PetscDSGetBdJacobian" 831 PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g, 832 void (**g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g0[]), 833 void (**g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g1[]), 834 void (**g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g2[]), 835 void (**g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g3[])) 836 { 837 PetscFunctionBegin; 838 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 839 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); 840 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); 841 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];} 842 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];} 843 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];} 844 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];} 845 PetscFunctionReturn(0); 846 } 847 848 #undef __FUNCT__ 849 #define __FUNCT__ "PetscDSSetBdJacobian" 850 PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g, 851 void (*g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g0[]), 852 void (*g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g1[]), 853 void (*g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g2[]), 854 void (*g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g3[])) 855 { 856 PetscErrorCode ierr; 857 858 PetscFunctionBegin; 859 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 860 if (g0) PetscValidFunction(g0, 4); 861 if (g1) PetscValidFunction(g1, 5); 862 if (g2) PetscValidFunction(g2, 6); 863 if (g3) PetscValidFunction(g3, 7); 864 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 865 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 866 ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr); 867 prob->gBd[(f*prob->Nf + g)*4+0] = g0; 868 prob->gBd[(f*prob->Nf + g)*4+1] = g1; 869 prob->gBd[(f*prob->Nf + g)*4+2] = g2; 870 prob->gBd[(f*prob->Nf + g)*4+3] = g3; 871 PetscFunctionReturn(0); 872 } 873 874 #undef __FUNCT__ 875 #define __FUNCT__ "PetscDSGetFieldOffset" 876 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off) 877 { 878 PetscInt g; 879 PetscErrorCode ierr; 880 881 PetscFunctionBegin; 882 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 883 PetscValidPointer(off, 3); 884 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); 885 *off = 0; 886 for (g = 0; g < f; ++g) { 887 PetscFE fe = (PetscFE) prob->disc[g]; 888 PetscInt Nb, Nc; 889 890 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 891 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 892 *off += Nb*Nc; 893 } 894 PetscFunctionReturn(0); 895 } 896 897 #undef __FUNCT__ 898 #define __FUNCT__ "PetscDSGetBdFieldOffset" 899 PetscErrorCode PetscDSGetBdFieldOffset(PetscDS prob, PetscInt f, PetscInt *off) 900 { 901 PetscInt g; 902 PetscErrorCode ierr; 903 904 PetscFunctionBegin; 905 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 906 PetscValidPointer(off, 3); 907 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); 908 *off = 0; 909 for (g = 0; g < f; ++g) { 910 PetscFE fe = (PetscFE) prob->discBd[g]; 911 PetscInt Nb, Nc; 912 913 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 914 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 915 *off += Nb*Nc; 916 } 917 PetscFunctionReturn(0); 918 } 919 920 #undef __FUNCT__ 921 #define __FUNCT__ "PetscDSGetTabulation" 922 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer) 923 { 924 PetscErrorCode ierr; 925 926 PetscFunctionBegin; 927 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 928 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 929 if (basis) {PetscValidPointer(basis, 2); *basis = prob->basis;} 930 if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;} 931 PetscFunctionReturn(0); 932 } 933 934 #undef __FUNCT__ 935 #define __FUNCT__ "PetscDSGetBdTabulation" 936 PetscErrorCode PetscDSGetBdTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer) 937 { 938 PetscErrorCode ierr; 939 940 PetscFunctionBegin; 941 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 942 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 943 if (basis) {PetscValidPointer(basis, 2); *basis = prob->basisBd;} 944 if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerBd;} 945 PetscFunctionReturn(0); 946 } 947 948 #undef __FUNCT__ 949 #define __FUNCT__ "PetscDSGetEvaluationArrays" 950 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x) 951 { 952 PetscErrorCode ierr; 953 954 PetscFunctionBegin; 955 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 956 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 957 if (u) {PetscValidPointer(u, 2); *u = prob->u;} 958 if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;} 959 if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;} 960 PetscFunctionReturn(0); 961 } 962 963 #undef __FUNCT__ 964 #define __FUNCT__ "PetscDSGetWeakFormArrays" 965 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3) 966 { 967 PetscErrorCode ierr; 968 969 PetscFunctionBegin; 970 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 971 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 972 if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;} 973 if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;} 974 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;} 975 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;} 976 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;} 977 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;} 978 PetscFunctionReturn(0); 979 } 980 981 #undef __FUNCT__ 982 #define __FUNCT__ "PetscDSGetRefCoordArrays" 983 PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer) 984 { 985 PetscErrorCode ierr; 986 987 PetscFunctionBegin; 988 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 989 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 990 if (x) {PetscValidPointer(x, 2); *x = prob->x;} 991 if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;} 992 PetscFunctionReturn(0); 993 } 994 995 #undef __FUNCT__ 996 #define __FUNCT__ "PetscDSDestroy_Basic" 997 PetscErrorCode PetscDSDestroy_Basic(PetscDS prob) 998 { 999 PetscFunctionBegin; 1000 PetscFunctionReturn(0); 1001 } 1002 1003 #undef __FUNCT__ 1004 #define __FUNCT__ "PetscDSInitialize_Basic" 1005 PetscErrorCode PetscDSInitialize_Basic(PetscDS prob) 1006 { 1007 PetscFunctionBegin; 1008 prob->ops->setfromoptions = NULL; 1009 prob->ops->setup = NULL; 1010 prob->ops->view = NULL; 1011 prob->ops->destroy = PetscDSDestroy_Basic; 1012 PetscFunctionReturn(0); 1013 } 1014 1015 /*MC 1016 PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions 1017 1018 Level: intermediate 1019 1020 .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType() 1021 M*/ 1022 1023 #undef __FUNCT__ 1024 #define __FUNCT__ "PetscDSCreate_Basic" 1025 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob) 1026 { 1027 PetscDS_Basic *b; 1028 PetscErrorCode ierr; 1029 1030 PetscFunctionBegin; 1031 PetscValidHeaderSpecific(prob, PETSCSPACE_CLASSID, 1); 1032 ierr = PetscNewLog(prob, &b);CHKERRQ(ierr); 1033 prob->data = b; 1034 1035 ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr); 1036 PetscFunctionReturn(0); 1037 } 1038