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