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 .seealso: PetscDSRegisterAll(), PetscDSRegisterDestroy() 53 54 @*/ 55 PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS)) 56 { 57 PetscErrorCode ierr; 58 59 PetscFunctionBegin; 60 ierr = PetscFunctionListAdd(&PetscDSList, sname, function);CHKERRQ(ierr); 61 PetscFunctionReturn(0); 62 } 63 64 /*@C 65 PetscDSSetType - Builds a particular PetscDS 66 67 Collective on prob 68 69 Input Parameters: 70 + prob - The PetscDS object 71 - name - The kind of system 72 73 Options Database Key: 74 . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types 75 76 Level: intermediate 77 78 Not available from Fortran 79 80 .seealso: PetscDSGetType(), PetscDSCreate() 81 @*/ 82 PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name) 83 { 84 PetscErrorCode (*r)(PetscDS); 85 PetscBool match; 86 PetscErrorCode ierr; 87 88 PetscFunctionBegin; 89 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 90 ierr = PetscObjectTypeCompare((PetscObject) prob, name, &match);CHKERRQ(ierr); 91 if (match) PetscFunctionReturn(0); 92 93 ierr = PetscDSRegisterAll();CHKERRQ(ierr); 94 ierr = PetscFunctionListFind(PetscDSList, name, &r);CHKERRQ(ierr); 95 if (!r) SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name); 96 97 if (prob->ops->destroy) { 98 ierr = (*prob->ops->destroy)(prob);CHKERRQ(ierr); 99 prob->ops->destroy = NULL; 100 } 101 ierr = (*r)(prob);CHKERRQ(ierr); 102 ierr = PetscObjectChangeTypeName((PetscObject) prob, name);CHKERRQ(ierr); 103 PetscFunctionReturn(0); 104 } 105 106 /*@C 107 PetscDSGetType - Gets the PetscDS type name (as a string) from the object. 108 109 Not Collective 110 111 Input Parameter: 112 . prob - The PetscDS 113 114 Output Parameter: 115 . name - The PetscDS type name 116 117 Level: intermediate 118 119 Not available from Fortran 120 121 .seealso: PetscDSSetType(), PetscDSCreate() 122 @*/ 123 PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name) 124 { 125 PetscErrorCode ierr; 126 127 PetscFunctionBegin; 128 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 129 PetscValidPointer(name, 2); 130 ierr = PetscDSRegisterAll();CHKERRQ(ierr); 131 *name = ((PetscObject) prob)->type_name; 132 PetscFunctionReturn(0); 133 } 134 135 static PetscErrorCode PetscDSView_Ascii(PetscDS prob, PetscViewer viewer) 136 { 137 PetscViewerFormat format; 138 const PetscScalar *constants; 139 PetscInt numConstants, f; 140 PetscErrorCode ierr; 141 142 PetscFunctionBegin; 143 ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 144 ierr = PetscViewerASCIIPrintf(viewer, "Discrete System with %d fields\n", prob->Nf);CHKERRQ(ierr); 145 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 146 ierr = PetscViewerASCIIPrintf(viewer, " cell total dim %D total comp %D\n", prob->totDim, prob->totComp);CHKERRQ(ierr); 147 if (prob->isHybrid) {ierr = PetscViewerASCIIPrintf(viewer, " hybrid cell\n");CHKERRQ(ierr);} 148 for (f = 0; f < prob->Nf; ++f) { 149 DSBoundary b; 150 PetscObject obj; 151 PetscClassId id; 152 PetscQuadrature q; 153 const char *name; 154 PetscInt Nc, Nq, Nqc; 155 156 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 157 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 158 ierr = PetscObjectGetName(obj, &name);CHKERRQ(ierr); 159 ierr = PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>");CHKERRQ(ierr); 160 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);CHKERRQ(ierr); 161 if (id == PETSCFE_CLASSID) { 162 ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr); 163 ierr = PetscFEGetQuadrature((PetscFE) obj, &q);CHKERRQ(ierr); 164 ierr = PetscViewerASCIIPrintf(viewer, " FEM");CHKERRQ(ierr); 165 } else if (id == PETSCFV_CLASSID) { 166 ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr); 167 ierr = PetscFVGetQuadrature((PetscFV) obj, &q);CHKERRQ(ierr); 168 ierr = PetscViewerASCIIPrintf(viewer, " FVM");CHKERRQ(ierr); 169 } 170 else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f); 171 if (Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "%D components", Nc);CHKERRQ(ierr);} 172 else {ierr = PetscViewerASCIIPrintf(viewer, "%D component ", Nc);CHKERRQ(ierr);} 173 if (prob->implicit[f]) {ierr = PetscViewerASCIIPrintf(viewer, " (implicit)");CHKERRQ(ierr);} 174 else {ierr = PetscViewerASCIIPrintf(viewer, " (explicit)");CHKERRQ(ierr);} 175 if (q) { 176 ierr = PetscQuadratureGetData(q, NULL, &Nqc, &Nq, NULL, NULL);CHKERRQ(ierr); 177 ierr = PetscViewerASCIIPrintf(viewer, " (Nq %D Nqc %D)", Nq, Nqc);CHKERRQ(ierr); 178 } 179 ierr = PetscViewerASCIIPrintf(viewer, " %D-jet", prob->jetDegree[f]);CHKERRQ(ierr); 180 ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 181 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);CHKERRQ(ierr); 182 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 183 if (id == PETSCFE_CLASSID) {ierr = PetscFEView((PetscFE) obj, viewer);CHKERRQ(ierr);} 184 else if (id == PETSCFV_CLASSID) {ierr = PetscFVView((PetscFV) obj, viewer);CHKERRQ(ierr);} 185 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 186 187 for (b = prob->boundary; b; b = b->next) { 188 PetscInt c, i; 189 190 if (b->field != f) continue; 191 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 192 ierr = PetscViewerASCIIPrintf(viewer, "Boundary %s (%s) %s\n", b->name, b->lname, DMBoundaryConditionTypes[b->type]);CHKERRQ(ierr); 193 if (!b->Nc) { 194 ierr = PetscViewerASCIIPrintf(viewer, " all components\n");CHKERRQ(ierr); 195 } else { 196 ierr = PetscViewerASCIIPrintf(viewer, " components: ");CHKERRQ(ierr); 197 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);CHKERRQ(ierr); 198 for (c = 0; c < b->Nc; ++c) { 199 if (c > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 200 ierr = PetscViewerASCIIPrintf(viewer, "%D", b->comps[c]);CHKERRQ(ierr); 201 } 202 ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 203 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);CHKERRQ(ierr); 204 } 205 ierr = PetscViewerASCIIPrintf(viewer, " values: ");CHKERRQ(ierr); 206 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);CHKERRQ(ierr); 207 for (i = 0; i < b->Nv; ++i) { 208 if (i > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 209 ierr = PetscViewerASCIIPrintf(viewer, "%D", b->values[i]);CHKERRQ(ierr); 210 } 211 ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 212 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);CHKERRQ(ierr); 213 ierr = PetscWeakFormView(b->wf, viewer);CHKERRQ(ierr); 214 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 215 } 216 } 217 ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 218 if (numConstants) { 219 ierr = PetscViewerASCIIPrintf(viewer, "%D constants\n", numConstants);CHKERRQ(ierr); 220 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 221 for (f = 0; f < numConstants; ++f) {ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double) PetscRealPart(constants[f]));CHKERRQ(ierr);} 222 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 223 } 224 ierr = PetscWeakFormView(prob->wf, viewer);CHKERRQ(ierr); 225 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 226 PetscFunctionReturn(0); 227 } 228 229 /*@C 230 PetscDSViewFromOptions - View from Options 231 232 Collective on PetscDS 233 234 Input Parameters: 235 + A - the PetscDS object 236 . obj - Optional object 237 - name - command line option 238 239 Level: intermediate 240 .seealso: PetscDS, PetscDSView, PetscObjectViewFromOptions(), PetscDSCreate() 241 @*/ 242 PetscErrorCode PetscDSViewFromOptions(PetscDS A,PetscObject obj,const char name[]) 243 { 244 PetscErrorCode ierr; 245 246 PetscFunctionBegin; 247 PetscValidHeaderSpecific(A,PETSCDS_CLASSID,1); 248 ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr); 249 PetscFunctionReturn(0); 250 } 251 252 /*@C 253 PetscDSView - Views a PetscDS 254 255 Collective on prob 256 257 Input Parameter: 258 + prob - the PetscDS object to view 259 - v - the viewer 260 261 Level: developer 262 263 .seealso PetscDSDestroy() 264 @*/ 265 PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v) 266 { 267 PetscBool iascii; 268 PetscErrorCode ierr; 269 270 PetscFunctionBegin; 271 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 272 if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);} 273 else {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);} 274 ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 275 if (iascii) {ierr = PetscDSView_Ascii(prob, v);CHKERRQ(ierr);} 276 if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);} 277 PetscFunctionReturn(0); 278 } 279 280 /*@ 281 PetscDSSetFromOptions - sets parameters in a PetscDS from the options database 282 283 Collective on prob 284 285 Input Parameter: 286 . prob - the PetscDS object to set options for 287 288 Options Database: 289 + -petscds_type <type> : Set the DS type 290 . -petscds_view <view opt> : View the DS 291 . -petscds_jac_pre : Turn formation of a separate Jacobian preconditioner on and off 292 . -bc_<name> <ids> : Specify a list of label ids for a boundary condition 293 - -bc_<name>_comp <comps> : Specify a list of field components to constrain for a boundary condition 294 295 Level: developer 296 297 .seealso PetscDSView() 298 @*/ 299 PetscErrorCode PetscDSSetFromOptions(PetscDS prob) 300 { 301 DSBoundary b; 302 const char *defaultType; 303 char name[256]; 304 PetscBool flg; 305 PetscErrorCode ierr; 306 307 PetscFunctionBegin; 308 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 309 if (!((PetscObject) prob)->type_name) { 310 defaultType = PETSCDSBASIC; 311 } else { 312 defaultType = ((PetscObject) prob)->type_name; 313 } 314 ierr = PetscDSRegisterAll();CHKERRQ(ierr); 315 316 ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr); 317 for (b = prob->boundary; b; b = b->next) { 318 char optname[1024]; 319 PetscInt ids[1024], len = 1024; 320 PetscBool flg; 321 322 ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr); 323 ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr); 324 ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr); 325 if (flg) { 326 b->Nv = len; 327 ierr = PetscFree(b->values);CHKERRQ(ierr); 328 ierr = PetscMalloc1(len, &b->values);CHKERRQ(ierr); 329 ierr = PetscArraycpy(b->values, ids, len);CHKERRQ(ierr); 330 ierr = PetscWeakFormRewriteKeys(b->wf, b->label, len, b->values);CHKERRQ(ierr); 331 } 332 len = 1024; 333 ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s_comp", b->name);CHKERRQ(ierr); 334 ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr); 335 ierr = PetscOptionsIntArray(optname, "List of boundary field components", "", ids, &len, &flg);CHKERRQ(ierr); 336 if (flg) { 337 b->Nc = len; 338 ierr = PetscFree(b->comps);CHKERRQ(ierr); 339 ierr = PetscMalloc1(len, &b->comps);CHKERRQ(ierr); 340 ierr = PetscArraycpy(b->comps, ids, len);CHKERRQ(ierr); 341 } 342 } 343 ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr); 344 if (flg) { 345 ierr = PetscDSSetType(prob, name);CHKERRQ(ierr); 346 } else if (!((PetscObject) prob)->type_name) { 347 ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr); 348 } 349 ierr = PetscOptionsBool("-petscds_jac_pre", "Discrete System", "PetscDSUseJacobianPreconditioner", prob->useJacPre, &prob->useJacPre, &flg);CHKERRQ(ierr); 350 if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);} 351 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 352 ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) prob);CHKERRQ(ierr); 353 ierr = PetscOptionsEnd();CHKERRQ(ierr); 354 if (prob->Nf) {ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr);} 355 PetscFunctionReturn(0); 356 } 357 358 /*@C 359 PetscDSSetUp - Construct data structures for the PetscDS 360 361 Collective on prob 362 363 Input Parameter: 364 . prob - the PetscDS object to setup 365 366 Level: developer 367 368 .seealso PetscDSView(), PetscDSDestroy() 369 @*/ 370 PetscErrorCode PetscDSSetUp(PetscDS prob) 371 { 372 const PetscInt Nf = prob->Nf; 373 PetscBool hasH = PETSC_FALSE; 374 PetscInt dim, dimEmbed, NbMax = 0, NcMax = 0, NqMax = 0, NsMax = 1, f; 375 PetscErrorCode ierr; 376 377 PetscFunctionBegin; 378 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 379 if (prob->setup) PetscFunctionReturn(0); 380 /* Calculate sizes */ 381 ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr); 382 ierr = PetscDSGetCoordinateDimension(prob, &dimEmbed);CHKERRQ(ierr); 383 prob->totDim = prob->totComp = 0; 384 ierr = PetscMalloc2(Nf,&prob->Nc,Nf,&prob->Nb);CHKERRQ(ierr); 385 ierr = PetscCalloc2(Nf+1,&prob->off,Nf+1,&prob->offDer);CHKERRQ(ierr); 386 ierr = PetscMalloc2(Nf,&prob->T,Nf,&prob->Tf);CHKERRQ(ierr); 387 for (f = 0; f < Nf; ++f) { 388 PetscObject obj; 389 PetscClassId id; 390 PetscQuadrature q = NULL; 391 PetscInt Nq = 0, Nb, Nc; 392 393 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 394 if (prob->jetDegree[f] > 1) hasH = PETSC_TRUE; 395 if (!obj) { 396 /* Empty mesh */ 397 Nb = Nc = 0; 398 prob->T[f] = prob->Tf[f] = NULL; 399 } else { 400 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 401 if (id == PETSCFE_CLASSID) { 402 PetscFE fe = (PetscFE) obj; 403 404 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 405 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 406 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 407 ierr = PetscFEGetCellTabulation(fe, prob->jetDegree[f], &prob->T[f]);CHKERRQ(ierr); 408 ierr = PetscFEGetFaceTabulation(fe, prob->jetDegree[f], &prob->Tf[f]);CHKERRQ(ierr); 409 } else if (id == PETSCFV_CLASSID) { 410 PetscFV fv = (PetscFV) obj; 411 412 ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr); 413 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 414 Nb = Nc; 415 ierr = PetscFVGetCellTabulation(fv, &prob->T[f]);CHKERRQ(ierr); 416 /* TODO: should PetscFV also have face tabulation? Otherwise there will be a null pointer in prob->basisFace */ 417 } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 418 } 419 prob->Nc[f] = Nc; 420 prob->Nb[f] = Nb; 421 prob->off[f+1] = Nc + prob->off[f]; 422 prob->offDer[f+1] = Nc*dim + prob->offDer[f]; 423 if (q) {ierr = PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);} 424 NqMax = PetscMax(NqMax, Nq); 425 NbMax = PetscMax(NbMax, Nb); 426 NcMax = PetscMax(NcMax, Nc); 427 prob->totDim += Nb; 428 prob->totComp += Nc; 429 /* There are two faces for all fields but the cohesive field on a hybrid cell */ 430 if (prob->isHybrid && (f < Nf-1)) prob->totDim += Nb; 431 } 432 /* Allocate works space */ 433 if (prob->isHybrid) NsMax = 2; 434 ierr = PetscMalloc3(NsMax*prob->totComp,&prob->u,NsMax*prob->totComp,&prob->u_t,NsMax*prob->totComp*dimEmbed + (hasH ? NsMax*prob->totComp*dimEmbed*dimEmbed : 0),&prob->u_x);CHKERRQ(ierr); 435 ierr = PetscMalloc5(dimEmbed,&prob->x,NbMax*NcMax,&prob->basisReal,NbMax*NcMax*dimEmbed,&prob->basisDerReal,NbMax*NcMax,&prob->testReal,NbMax*NcMax*dimEmbed,&prob->testDerReal);CHKERRQ(ierr); 436 ierr = PetscMalloc6(NsMax*NqMax*NcMax,&prob->f0,NsMax*NqMax*NcMax*dimEmbed,&prob->f1, 437 NsMax*NsMax*NqMax*NcMax*NcMax,&prob->g0,NsMax*NsMax*NqMax*NcMax*NcMax*dimEmbed,&prob->g1, 438 NsMax*NsMax*NqMax*NcMax*NcMax*dimEmbed,&prob->g2,NsMax*NsMax*NqMax*NcMax*NcMax*dimEmbed*dimEmbed,&prob->g3);CHKERRQ(ierr); 439 if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);} 440 prob->setup = PETSC_TRUE; 441 PetscFunctionReturn(0); 442 } 443 444 static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob) 445 { 446 PetscErrorCode ierr; 447 448 PetscFunctionBegin; 449 ierr = PetscFree2(prob->Nc,prob->Nb);CHKERRQ(ierr); 450 ierr = PetscFree2(prob->off,prob->offDer);CHKERRQ(ierr); 451 ierr = PetscFree2(prob->T,prob->Tf);CHKERRQ(ierr); 452 ierr = PetscFree3(prob->u,prob->u_t,prob->u_x);CHKERRQ(ierr); 453 ierr = PetscFree5(prob->x,prob->basisReal, prob->basisDerReal,prob->testReal,prob->testDerReal);CHKERRQ(ierr); 454 ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr); 455 PetscFunctionReturn(0); 456 } 457 458 static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew) 459 { 460 PetscObject *tmpd; 461 PetscBool *tmpi; 462 PetscInt *tmpk; 463 PetscPointFunc *tmpup; 464 PetscSimplePointFunc *tmpexactSol, *tmpexactSol_t; 465 void **tmpexactCtx, **tmpexactCtx_t; 466 void **tmpctx; 467 PetscInt Nf = prob->Nf, f; 468 PetscErrorCode ierr; 469 470 PetscFunctionBegin; 471 if (Nf >= NfNew) PetscFunctionReturn(0); 472 prob->setup = PETSC_FALSE; 473 ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr); 474 ierr = PetscMalloc3(NfNew, &tmpd, NfNew, &tmpi, NfNew, &tmpk);CHKERRQ(ierr); 475 for (f = 0; f < Nf; ++f) {tmpd[f] = prob->disc[f]; tmpi[f] = prob->implicit[f]; tmpk[f] = prob->jetDegree[f];} 476 for (f = Nf; f < NfNew; ++f) {tmpd[f] = NULL; tmpi[f] = PETSC_TRUE; tmpk[f] = 1;} 477 ierr = PetscFree3(prob->disc, prob->implicit, prob->jetDegree);CHKERRQ(ierr); 478 ierr = PetscWeakFormSetNumFields(prob->wf, NfNew);CHKERRQ(ierr); 479 prob->Nf = NfNew; 480 prob->disc = tmpd; 481 prob->implicit = tmpi; 482 prob->jetDegree = tmpk; 483 ierr = PetscCalloc2(NfNew, &tmpup, NfNew, &tmpctx);CHKERRQ(ierr); 484 for (f = 0; f < Nf; ++f) tmpup[f] = prob->update[f]; 485 for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f]; 486 for (f = Nf; f < NfNew; ++f) tmpup[f] = NULL; 487 for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL; 488 ierr = PetscFree2(prob->update, prob->ctx);CHKERRQ(ierr); 489 prob->update = tmpup; 490 prob->ctx = tmpctx; 491 ierr = PetscCalloc4(NfNew, &tmpexactSol, NfNew, &tmpexactCtx, NfNew, &tmpexactSol_t, NfNew, &tmpexactCtx_t);CHKERRQ(ierr); 492 for (f = 0; f < Nf; ++f) tmpexactSol[f] = prob->exactSol[f]; 493 for (f = 0; f < Nf; ++f) tmpexactCtx[f] = prob->exactCtx[f]; 494 for (f = 0; f < Nf; ++f) tmpexactSol_t[f] = prob->exactSol_t[f]; 495 for (f = 0; f < Nf; ++f) tmpexactCtx_t[f] = prob->exactCtx_t[f]; 496 for (f = Nf; f < NfNew; ++f) tmpexactSol[f] = NULL; 497 for (f = Nf; f < NfNew; ++f) tmpexactCtx[f] = NULL; 498 for (f = Nf; f < NfNew; ++f) tmpexactSol_t[f] = NULL; 499 for (f = Nf; f < NfNew; ++f) tmpexactCtx_t[f] = NULL; 500 ierr = PetscFree4(prob->exactSol, prob->exactCtx, prob->exactSol_t, prob->exactCtx_t);CHKERRQ(ierr); 501 prob->exactSol = tmpexactSol; 502 prob->exactCtx = tmpexactCtx; 503 prob->exactSol_t = tmpexactSol_t; 504 prob->exactCtx_t = tmpexactCtx_t; 505 PetscFunctionReturn(0); 506 } 507 508 /*@ 509 PetscDSDestroy - Destroys a PetscDS object 510 511 Collective on prob 512 513 Input Parameter: 514 . prob - the PetscDS object to destroy 515 516 Level: developer 517 518 .seealso PetscDSView() 519 @*/ 520 PetscErrorCode PetscDSDestroy(PetscDS *ds) 521 { 522 PetscInt f; 523 PetscErrorCode ierr; 524 525 PetscFunctionBegin; 526 if (!*ds) PetscFunctionReturn(0); 527 PetscValidHeaderSpecific((*ds), PETSCDS_CLASSID, 1); 528 529 if (--((PetscObject)(*ds))->refct > 0) {*ds = NULL; PetscFunctionReturn(0);} 530 ((PetscObject) (*ds))->refct = 0; 531 if ((*ds)->subprobs) { 532 PetscInt dim, d; 533 534 ierr = PetscDSGetSpatialDimension(*ds, &dim);CHKERRQ(ierr); 535 for (d = 0; d < dim; ++d) {ierr = PetscDSDestroy(&(*ds)->subprobs[d]);CHKERRQ(ierr);} 536 } 537 ierr = PetscFree((*ds)->subprobs);CHKERRQ(ierr); 538 ierr = PetscDSDestroyStructs_Static(*ds);CHKERRQ(ierr); 539 for (f = 0; f < (*ds)->Nf; ++f) { 540 ierr = PetscObjectDereference((*ds)->disc[f]);CHKERRQ(ierr); 541 } 542 ierr = PetscFree3((*ds)->disc, (*ds)->implicit, (*ds)->jetDegree);CHKERRQ(ierr); 543 ierr = PetscWeakFormDestroy(&(*ds)->wf);CHKERRQ(ierr); 544 ierr = PetscFree2((*ds)->update,(*ds)->ctx);CHKERRQ(ierr); 545 ierr = PetscFree4((*ds)->exactSol,(*ds)->exactCtx,(*ds)->exactSol_t,(*ds)->exactCtx_t);CHKERRQ(ierr); 546 if ((*ds)->ops->destroy) {ierr = (*(*ds)->ops->destroy)(*ds);CHKERRQ(ierr);} 547 ierr = PetscDSDestroyBoundary(*ds);CHKERRQ(ierr); 548 ierr = PetscFree((*ds)->constants);CHKERRQ(ierr); 549 ierr = PetscHeaderDestroy(ds);CHKERRQ(ierr); 550 PetscFunctionReturn(0); 551 } 552 553 /*@ 554 PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType(). 555 556 Collective 557 558 Input Parameter: 559 . comm - The communicator for the PetscDS object 560 561 Output Parameter: 562 . ds - The PetscDS object 563 564 Level: beginner 565 566 .seealso: PetscDSSetType(), PETSCDSBASIC 567 @*/ 568 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *ds) 569 { 570 PetscDS p; 571 PetscErrorCode ierr; 572 573 PetscFunctionBegin; 574 PetscValidPointer(ds, 2); 575 *ds = NULL; 576 ierr = PetscDSInitializePackage();CHKERRQ(ierr); 577 578 ierr = PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr); 579 580 p->Nf = 0; 581 p->setup = PETSC_FALSE; 582 p->numConstants = 0; 583 p->constants = NULL; 584 p->dimEmbed = -1; 585 p->useJacPre = PETSC_TRUE; 586 ierr = PetscWeakFormCreate(comm, &p->wf);CHKERRQ(ierr); 587 588 *ds = p; 589 PetscFunctionReturn(0); 590 } 591 592 /*@ 593 PetscDSGetNumFields - Returns the number of fields in the DS 594 595 Not collective 596 597 Input Parameter: 598 . prob - The PetscDS object 599 600 Output Parameter: 601 . Nf - The number of fields 602 603 Level: beginner 604 605 .seealso: PetscDSGetSpatialDimension(), PetscDSCreate() 606 @*/ 607 PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf) 608 { 609 PetscFunctionBegin; 610 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 611 PetscValidPointer(Nf, 2); 612 *Nf = prob->Nf; 613 PetscFunctionReturn(0); 614 } 615 616 /*@ 617 PetscDSGetSpatialDimension - Returns the spatial dimension of the DS, meaning the topological dimension of the discretizations 618 619 Not collective 620 621 Input Parameter: 622 . prob - The PetscDS object 623 624 Output Parameter: 625 . dim - The spatial dimension 626 627 Level: beginner 628 629 .seealso: PetscDSGetCoordinateDimension(), PetscDSGetNumFields(), PetscDSCreate() 630 @*/ 631 PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim) 632 { 633 PetscErrorCode ierr; 634 635 PetscFunctionBegin; 636 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 637 PetscValidPointer(dim, 2); 638 *dim = 0; 639 if (prob->Nf) { 640 PetscObject obj; 641 PetscClassId id; 642 643 ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr); 644 if (obj) { 645 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 646 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);} 647 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);} 648 else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0); 649 } 650 } 651 PetscFunctionReturn(0); 652 } 653 654 /*@ 655 PetscDSGetCoordinateDimension - Returns the coordinate dimension of the DS, meaning the dimension of the space into which the discretiaztions are embedded 656 657 Not collective 658 659 Input Parameter: 660 . prob - The PetscDS object 661 662 Output Parameter: 663 . dimEmbed - The coordinate dimension 664 665 Level: beginner 666 667 .seealso: PetscDSSetCoordinateDimension(), PetscDSGetSpatialDimension(), PetscDSGetNumFields(), PetscDSCreate() 668 @*/ 669 PetscErrorCode PetscDSGetCoordinateDimension(PetscDS prob, PetscInt *dimEmbed) 670 { 671 PetscFunctionBegin; 672 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 673 PetscValidPointer(dimEmbed, 2); 674 if (prob->dimEmbed < 0) SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONGSTATE, "No coordinate dimension set for this DS"); 675 *dimEmbed = prob->dimEmbed; 676 PetscFunctionReturn(0); 677 } 678 679 /*@ 680 PetscDSSetCoordinateDimension - Set the coordinate dimension of the DS, meaning the dimension of the space into which the discretiaztions are embedded 681 682 Logically collective on prob 683 684 Input Parameters: 685 + prob - The PetscDS object 686 - dimEmbed - The coordinate dimension 687 688 Level: beginner 689 690 .seealso: PetscDSGetCoordinateDimension(), PetscDSGetSpatialDimension(), PetscDSGetNumFields(), PetscDSCreate() 691 @*/ 692 PetscErrorCode PetscDSSetCoordinateDimension(PetscDS prob, PetscInt dimEmbed) 693 { 694 PetscFunctionBegin; 695 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 696 if (dimEmbed < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coordinate dimension must be non-negative, not %D", dimEmbed); 697 prob->dimEmbed = dimEmbed; 698 PetscFunctionReturn(0); 699 } 700 701 /*@ 702 PetscDSGetHybrid - Returns the flag for a hybrid (cohesive) cell 703 704 Not collective 705 706 Input Parameter: 707 . prob - The PetscDS object 708 709 Output Parameter: 710 . isHybrid - The flag 711 712 Level: developer 713 714 .seealso: PetscDSSetHybrid(), PetscDSCreate() 715 @*/ 716 PetscErrorCode PetscDSGetHybrid(PetscDS prob, PetscBool *isHybrid) 717 { 718 PetscFunctionBegin; 719 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 720 PetscValidPointer(isHybrid, 2); 721 *isHybrid = prob->isHybrid; 722 PetscFunctionReturn(0); 723 } 724 725 /*@ 726 PetscDSSetHybrid - Set the flag for a hybrid (cohesive) cell 727 728 Not collective 729 730 Input Parameters: 731 + prob - The PetscDS object 732 - isHybrid - The flag 733 734 Level: developer 735 736 .seealso: PetscDSGetHybrid(), PetscDSCreate() 737 @*/ 738 PetscErrorCode PetscDSSetHybrid(PetscDS prob, PetscBool isHybrid) 739 { 740 PetscFunctionBegin; 741 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 742 prob->isHybrid = isHybrid; 743 PetscFunctionReturn(0); 744 } 745 746 /*@ 747 PetscDSGetTotalDimension - Returns the total size of the approximation space for this system 748 749 Not collective 750 751 Input Parameter: 752 . prob - The PetscDS object 753 754 Output Parameter: 755 . dim - The total problem dimension 756 757 Level: beginner 758 759 .seealso: PetscDSGetNumFields(), PetscDSCreate() 760 @*/ 761 PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim) 762 { 763 PetscErrorCode ierr; 764 765 PetscFunctionBegin; 766 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 767 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 768 PetscValidPointer(dim, 2); 769 *dim = prob->totDim; 770 PetscFunctionReturn(0); 771 } 772 773 /*@ 774 PetscDSGetTotalComponents - Returns the total number of components in this system 775 776 Not collective 777 778 Input Parameter: 779 . prob - The PetscDS object 780 781 Output Parameter: 782 . dim - The total number of components 783 784 Level: beginner 785 786 .seealso: PetscDSGetNumFields(), PetscDSCreate() 787 @*/ 788 PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc) 789 { 790 PetscErrorCode ierr; 791 792 PetscFunctionBegin; 793 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 794 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 795 PetscValidPointer(Nc, 2); 796 *Nc = prob->totComp; 797 PetscFunctionReturn(0); 798 } 799 800 /*@ 801 PetscDSGetDiscretization - Returns the discretization object for the given field 802 803 Not collective 804 805 Input Parameters: 806 + prob - The PetscDS object 807 - f - The field number 808 809 Output Parameter: 810 . disc - The discretization object 811 812 Level: beginner 813 814 .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 815 @*/ 816 PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc) 817 { 818 PetscFunctionBeginHot; 819 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 820 PetscValidPointer(disc, 3); 821 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); 822 *disc = prob->disc[f]; 823 PetscFunctionReturn(0); 824 } 825 826 /*@ 827 PetscDSSetDiscretization - Sets the discretization object for the given field 828 829 Not collective 830 831 Input Parameters: 832 + prob - The PetscDS object 833 . f - The field number 834 - disc - The discretization object 835 836 Level: beginner 837 838 .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 839 @*/ 840 PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc) 841 { 842 PetscErrorCode ierr; 843 844 PetscFunctionBegin; 845 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 846 if (disc) PetscValidPointer(disc, 3); 847 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 848 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 849 ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr); 850 prob->disc[f] = disc; 851 ierr = PetscObjectReference(disc);CHKERRQ(ierr); 852 if (disc) { 853 PetscClassId id; 854 855 ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr); 856 if (id == PETSCFE_CLASSID) { 857 ierr = PetscDSSetImplicit(prob, f, PETSC_TRUE);CHKERRQ(ierr); 858 } else if (id == PETSCFV_CLASSID) { 859 ierr = PetscDSSetImplicit(prob, f, PETSC_FALSE);CHKERRQ(ierr); 860 } 861 ierr = PetscDSSetJetDegree(prob, f, 1);CHKERRQ(ierr); 862 } 863 PetscFunctionReturn(0); 864 } 865 866 /*@ 867 PetscDSGetWeakForm - Returns the weak form object 868 869 Not collective 870 871 Input Parameter: 872 . ds - The PetscDS object 873 874 Output Parameter: 875 . wf - The weak form object 876 877 Level: beginner 878 879 .seealso: PetscDSSetWeakForm(), PetscDSGetNumFields(), PetscDSCreate() 880 @*/ 881 PetscErrorCode PetscDSGetWeakForm(PetscDS ds, PetscWeakForm *wf) 882 { 883 PetscFunctionBegin; 884 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 885 PetscValidPointer(wf, 2); 886 *wf = ds->wf; 887 PetscFunctionReturn(0); 888 } 889 890 /*@ 891 PetscDSSetWeakForm - Sets the weak form object 892 893 Not collective 894 895 Input Parameters: 896 + ds - The PetscDS object 897 - wf - The weak form object 898 899 Level: beginner 900 901 .seealso: PetscDSGetWeakForm(), PetscDSGetNumFields(), PetscDSCreate() 902 @*/ 903 PetscErrorCode PetscDSSetWeakForm(PetscDS ds, PetscWeakForm wf) 904 { 905 PetscErrorCode ierr; 906 907 PetscFunctionBegin; 908 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 909 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 2); 910 ierr = PetscObjectDereference((PetscObject) ds->wf);CHKERRQ(ierr); 911 ds->wf = wf; 912 ierr = PetscObjectReference((PetscObject) wf);CHKERRQ(ierr); 913 ierr = PetscWeakFormSetNumFields(wf, ds->Nf);CHKERRQ(ierr); 914 PetscFunctionReturn(0); 915 } 916 917 /*@ 918 PetscDSAddDiscretization - Adds a discretization object 919 920 Not collective 921 922 Input Parameters: 923 + prob - The PetscDS object 924 - disc - The boundary discretization object 925 926 Level: beginner 927 928 .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 929 @*/ 930 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc) 931 { 932 PetscErrorCode ierr; 933 934 PetscFunctionBegin; 935 ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr); 936 PetscFunctionReturn(0); 937 } 938 939 /*@ 940 PetscDSGetQuadrature - Returns the quadrature, which must agree for all fields in the DS 941 942 Not collective 943 944 Input Parameter: 945 . prob - The PetscDS object 946 947 Output Parameter: 948 . q - The quadrature object 949 950 Level: intermediate 951 952 .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 953 @*/ 954 PetscErrorCode PetscDSGetQuadrature(PetscDS prob, PetscQuadrature *q) 955 { 956 PetscObject obj; 957 PetscClassId id; 958 PetscErrorCode ierr; 959 960 PetscFunctionBegin; 961 *q = NULL; 962 if (!prob->Nf) PetscFunctionReturn(0); 963 ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr); 964 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 965 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetQuadrature((PetscFE) obj, q);CHKERRQ(ierr);} 966 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetQuadrature((PetscFV) obj, q);CHKERRQ(ierr);} 967 else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0); 968 PetscFunctionReturn(0); 969 } 970 971 /*@ 972 PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX 973 974 Not collective 975 976 Input Parameters: 977 + prob - The PetscDS object 978 - f - The field number 979 980 Output Parameter: 981 . implicit - The flag indicating what kind of solve to use for this field 982 983 Level: developer 984 985 .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 986 @*/ 987 PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit) 988 { 989 PetscFunctionBegin; 990 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 991 PetscValidPointer(implicit, 3); 992 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); 993 *implicit = prob->implicit[f]; 994 PetscFunctionReturn(0); 995 } 996 997 /*@ 998 PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX 999 1000 Not collective 1001 1002 Input Parameters: 1003 + prob - The PetscDS object 1004 . f - The field number 1005 - implicit - The flag indicating what kind of solve to use for this field 1006 1007 Level: developer 1008 1009 .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 1010 @*/ 1011 PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit) 1012 { 1013 PetscFunctionBegin; 1014 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1015 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); 1016 prob->implicit[f] = implicit; 1017 PetscFunctionReturn(0); 1018 } 1019 1020 /*@ 1021 PetscDSGetJetDegree - Returns the highest derivative for this field equation, or the k-jet that the discretization needs to tabulate. 1022 1023 Not collective 1024 1025 Input Parameters: 1026 + ds - The PetscDS object 1027 - f - The field number 1028 1029 Output Parameter: 1030 . k - The highest derivative we need to tabulate 1031 1032 Level: developer 1033 1034 .seealso: PetscDSSetJetDegree(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 1035 @*/ 1036 PetscErrorCode PetscDSGetJetDegree(PetscDS ds, PetscInt f, PetscInt *k) 1037 { 1038 PetscFunctionBegin; 1039 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1040 PetscValidPointer(k, 3); 1041 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 1042 *k = ds->jetDegree[f]; 1043 PetscFunctionReturn(0); 1044 } 1045 1046 /*@ 1047 PetscDSSetJetDegree - Set the highest derivative for this field equation, or the k-jet that the discretization needs to tabulate. 1048 1049 Not collective 1050 1051 Input Parameters: 1052 + ds - The PetscDS object 1053 . f - The field number 1054 - k - The highest derivative we need to tabulate 1055 1056 Level: developer 1057 1058 .seealso: PetscDSGetJetDegree(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 1059 @*/ 1060 PetscErrorCode PetscDSSetJetDegree(PetscDS ds, PetscInt f, PetscInt k) 1061 { 1062 PetscFunctionBegin; 1063 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1064 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 1065 ds->jetDegree[f] = k; 1066 PetscFunctionReturn(0); 1067 } 1068 1069 PetscErrorCode PetscDSGetObjective(PetscDS ds, PetscInt f, 1070 void (**obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1071 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1072 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1073 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])) 1074 { 1075 PetscPointFunc *tmp; 1076 PetscInt n; 1077 PetscErrorCode ierr; 1078 1079 PetscFunctionBegin; 1080 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1081 PetscValidPointer(obj, 3); 1082 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 1083 ierr = PetscWeakFormGetObjective(ds->wf, NULL, 0, f, &n, &tmp);CHKERRQ(ierr); 1084 *obj = tmp ? tmp[0] : NULL; 1085 PetscFunctionReturn(0); 1086 } 1087 1088 PetscErrorCode PetscDSSetObjective(PetscDS ds, PetscInt f, 1089 void (*obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1090 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1091 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1092 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])) 1093 { 1094 PetscErrorCode ierr; 1095 1096 PetscFunctionBegin; 1097 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1098 if (obj) PetscValidFunction(obj, 3); 1099 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1100 ierr = PetscWeakFormSetIndexObjective(ds->wf, NULL, 0, f, 0, obj);CHKERRQ(ierr); 1101 PetscFunctionReturn(0); 1102 } 1103 1104 /*@C 1105 PetscDSGetResidual - Get the pointwise residual function for a given test field 1106 1107 Not collective 1108 1109 Input Parameters: 1110 + ds - The PetscDS 1111 - f - The test field number 1112 1113 Output Parameters: 1114 + f0 - integrand for the test function term 1115 - f1 - integrand for the test function gradient term 1116 1117 Note: We are using a first order FEM model for the weak form: 1118 1119 \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) 1120 1121 The calling sequence for the callbacks f0 and f1 is given by: 1122 1123 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1124 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1125 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1126 $ PetscReal t, const PetscReal x[], PetscScalar f0[]) 1127 1128 + dim - the spatial dimension 1129 . Nf - the number of fields 1130 . uOff - the offset into u[] and u_t[] for each field 1131 . uOff_x - the offset into u_x[] for each field 1132 . u - each field evaluated at the current point 1133 . u_t - the time derivative of each field evaluated at the current point 1134 . u_x - the gradient of each field evaluated at the current point 1135 . aOff - the offset into a[] and a_t[] for each auxiliary field 1136 . aOff_x - the offset into a_x[] for each auxiliary field 1137 . a - each auxiliary field evaluated at the current point 1138 . a_t - the time derivative of each auxiliary field evaluated at the current point 1139 . a_x - the gradient of auxiliary each field evaluated at the current point 1140 . t - current time 1141 . x - coordinates of the current point 1142 . numConstants - number of constant parameters 1143 . constants - constant parameters 1144 - f0 - output values at the current point 1145 1146 Level: intermediate 1147 1148 .seealso: PetscDSSetResidual() 1149 @*/ 1150 PetscErrorCode PetscDSGetResidual(PetscDS ds, PetscInt f, 1151 void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1152 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1153 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1154 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), 1155 void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1156 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1157 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1158 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 1159 { 1160 PetscPointFunc *tmp0, *tmp1; 1161 PetscInt n0, n1; 1162 PetscErrorCode ierr; 1163 1164 PetscFunctionBegin; 1165 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1166 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 1167 ierr = PetscWeakFormGetResidual(ds->wf, NULL, 0, f, &n0, &tmp0, &n1, &tmp1);CHKERRQ(ierr); 1168 *f0 = tmp0 ? tmp0[0] : NULL; 1169 *f1 = tmp1 ? tmp1[0] : NULL; 1170 PetscFunctionReturn(0); 1171 } 1172 1173 /*@C 1174 PetscDSSetResidual - Set the pointwise residual function for a given test field 1175 1176 Not collective 1177 1178 Input Parameters: 1179 + ds - The PetscDS 1180 . f - The test field number 1181 . f0 - integrand for the test function term 1182 - f1 - integrand for the test function gradient term 1183 1184 Note: We are using a first order FEM model for the weak form: 1185 1186 \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) 1187 1188 The calling sequence for the callbacks f0 and f1 is given by: 1189 1190 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1191 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1192 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1193 $ PetscReal t, const PetscReal x[], PetscScalar f0[]) 1194 1195 + dim - the spatial dimension 1196 . Nf - the number of fields 1197 . uOff - the offset into u[] and u_t[] for each field 1198 . uOff_x - the offset into u_x[] for each field 1199 . u - each field evaluated at the current point 1200 . u_t - the time derivative of each field evaluated at the current point 1201 . u_x - the gradient of each field evaluated at the current point 1202 . aOff - the offset into a[] and a_t[] for each auxiliary field 1203 . aOff_x - the offset into a_x[] for each auxiliary field 1204 . a - each auxiliary field evaluated at the current point 1205 . a_t - the time derivative of each auxiliary field evaluated at the current point 1206 . a_x - the gradient of auxiliary each field evaluated at the current point 1207 . t - current time 1208 . x - coordinates of the current point 1209 . numConstants - number of constant parameters 1210 . constants - constant parameters 1211 - f0 - output values at the current point 1212 1213 Level: intermediate 1214 1215 .seealso: PetscDSGetResidual() 1216 @*/ 1217 PetscErrorCode PetscDSSetResidual(PetscDS ds, PetscInt f, 1218 void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1219 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1220 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1221 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), 1222 void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1223 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1224 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1225 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 1226 { 1227 PetscErrorCode ierr; 1228 1229 PetscFunctionBegin; 1230 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1231 if (f0) PetscValidFunction(f0, 3); 1232 if (f1) PetscValidFunction(f1, 4); 1233 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1234 ierr = PetscWeakFormSetIndexResidual(ds->wf, NULL, 0, f, 0, f0, 0, f1);CHKERRQ(ierr); 1235 PetscFunctionReturn(0); 1236 } 1237 1238 /*@C 1239 PetscDSHasJacobian - Signals that Jacobian functions have been set 1240 1241 Not collective 1242 1243 Input Parameter: 1244 . prob - The PetscDS 1245 1246 Output Parameter: 1247 . hasJac - flag that pointwise function for the Jacobian has been set 1248 1249 Level: intermediate 1250 1251 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian() 1252 @*/ 1253 PetscErrorCode PetscDSHasJacobian(PetscDS ds, PetscBool *hasJac) 1254 { 1255 PetscErrorCode ierr; 1256 1257 PetscFunctionBegin; 1258 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1259 ierr = PetscWeakFormHasJacobian(ds->wf, hasJac);CHKERRQ(ierr); 1260 PetscFunctionReturn(0); 1261 } 1262 1263 /*@C 1264 PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field 1265 1266 Not collective 1267 1268 Input Parameters: 1269 + ds - The PetscDS 1270 . f - The test field number 1271 - g - The field number 1272 1273 Output Parameters: 1274 + g0 - integrand for the test and basis function term 1275 . g1 - integrand for the test function and basis function gradient term 1276 . g2 - integrand for the test function gradient and basis function term 1277 - g3 - integrand for the test function gradient and basis function gradient term 1278 1279 Note: We are using a first order FEM model for the weak form: 1280 1281 \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 1282 1283 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1284 1285 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1286 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1287 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1288 $ PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]) 1289 1290 + dim - the spatial dimension 1291 . Nf - the number of fields 1292 . uOff - the offset into u[] and u_t[] for each field 1293 . uOff_x - the offset into u_x[] for each field 1294 . u - each field evaluated at the current point 1295 . u_t - the time derivative of each field evaluated at the current point 1296 . u_x - the gradient of each field evaluated at the current point 1297 . aOff - the offset into a[] and a_t[] for each auxiliary field 1298 . aOff_x - the offset into a_x[] for each auxiliary field 1299 . a - each auxiliary field evaluated at the current point 1300 . a_t - the time derivative of each auxiliary field evaluated at the current point 1301 . a_x - the gradient of auxiliary each field evaluated at the current point 1302 . t - current time 1303 . u_tShift - the multiplier a for dF/dU_t 1304 . x - coordinates of the current point 1305 . numConstants - number of constant parameters 1306 . constants - constant parameters 1307 - g0 - output values at the current point 1308 1309 Level: intermediate 1310 1311 .seealso: PetscDSSetJacobian() 1312 @*/ 1313 PetscErrorCode PetscDSGetJacobian(PetscDS ds, PetscInt f, PetscInt g, 1314 void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1315 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1316 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1317 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1318 void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1319 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1320 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1321 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1322 void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1323 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1324 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1325 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1326 void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1327 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1328 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1329 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1330 { 1331 PetscPointJac *tmp0, *tmp1, *tmp2, *tmp3; 1332 PetscInt n0, n1, n2, n3; 1333 PetscErrorCode ierr; 1334 1335 PetscFunctionBegin; 1336 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1337 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 1338 if ((g < 0) || (g >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, ds->Nf); 1339 ierr = PetscWeakFormGetJacobian(ds->wf, NULL, 0, f, g, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3);CHKERRQ(ierr); 1340 *g0 = tmp0 ? tmp0[0] : NULL; 1341 *g1 = tmp1 ? tmp1[0] : NULL; 1342 *g2 = tmp2 ? tmp2[0] : NULL; 1343 *g3 = tmp3 ? tmp3[0] : NULL; 1344 PetscFunctionReturn(0); 1345 } 1346 1347 /*@C 1348 PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields 1349 1350 Not collective 1351 1352 Input Parameters: 1353 + ds - The PetscDS 1354 . f - The test field number 1355 . g - The field number 1356 . g0 - integrand for the test and basis function term 1357 . g1 - integrand for the test function and basis function gradient term 1358 . g2 - integrand for the test function gradient and basis function term 1359 - g3 - integrand for the test function gradient and basis function gradient term 1360 1361 Note: We are using a first order FEM model for the weak form: 1362 1363 \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 1364 1365 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1366 1367 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1368 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1369 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1370 $ PetscReal t, const PetscReal x[], PetscScalar g0[]) 1371 1372 + dim - the spatial dimension 1373 . Nf - the number of fields 1374 . uOff - the offset into u[] and u_t[] for each field 1375 . uOff_x - the offset into u_x[] for each field 1376 . u - each field evaluated at the current point 1377 . u_t - the time derivative of each field evaluated at the current point 1378 . u_x - the gradient of each field evaluated at the current point 1379 . aOff - the offset into a[] and a_t[] for each auxiliary field 1380 . aOff_x - the offset into a_x[] for each auxiliary field 1381 . a - each auxiliary field evaluated at the current point 1382 . a_t - the time derivative of each auxiliary field evaluated at the current point 1383 . a_x - the gradient of auxiliary each field evaluated at the current point 1384 . t - current time 1385 . u_tShift - the multiplier a for dF/dU_t 1386 . x - coordinates of the current point 1387 . numConstants - number of constant parameters 1388 . constants - constant parameters 1389 - g0 - output values at the current point 1390 1391 Level: intermediate 1392 1393 .seealso: PetscDSGetJacobian() 1394 @*/ 1395 PetscErrorCode PetscDSSetJacobian(PetscDS ds, PetscInt f, PetscInt g, 1396 void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1397 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1398 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1399 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1400 void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1401 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1402 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1403 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1404 void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1405 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1406 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1407 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1408 void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1409 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1410 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1411 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1412 { 1413 PetscErrorCode ierr; 1414 1415 PetscFunctionBegin; 1416 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1417 if (g0) PetscValidFunction(g0, 4); 1418 if (g1) PetscValidFunction(g1, 5); 1419 if (g2) PetscValidFunction(g2, 6); 1420 if (g3) PetscValidFunction(g3, 7); 1421 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1422 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 1423 ierr = PetscWeakFormSetIndexJacobian(ds->wf, NULL, 0, f, g, 0, g0, 0, g1, 0, g2, 0, g3);CHKERRQ(ierr); 1424 PetscFunctionReturn(0); 1425 } 1426 1427 /*@C 1428 PetscDSUseJacobianPreconditioner - Whether to construct a Jacobian preconditioner 1429 1430 Not collective 1431 1432 Input Parameters: 1433 + prob - The PetscDS 1434 - useJacPre - flag that enables construction of a Jacobian preconditioner 1435 1436 Level: intermediate 1437 1438 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian() 1439 @*/ 1440 PetscErrorCode PetscDSUseJacobianPreconditioner(PetscDS prob, PetscBool useJacPre) 1441 { 1442 PetscFunctionBegin; 1443 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1444 prob->useJacPre = useJacPre; 1445 PetscFunctionReturn(0); 1446 } 1447 1448 /*@C 1449 PetscDSHasJacobianPreconditioner - Signals that a Jacobian preconditioner matrix has been set 1450 1451 Not collective 1452 1453 Input Parameter: 1454 . prob - The PetscDS 1455 1456 Output Parameter: 1457 . hasJacPre - flag that pointwise function for Jacobian preconditioner matrix has been set 1458 1459 Level: intermediate 1460 1461 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian() 1462 @*/ 1463 PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS ds, PetscBool *hasJacPre) 1464 { 1465 PetscErrorCode ierr; 1466 1467 PetscFunctionBegin; 1468 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1469 *hasJacPre = PETSC_FALSE; 1470 if (!ds->useJacPre) PetscFunctionReturn(0); 1471 ierr = PetscWeakFormHasJacobianPreconditioner(ds->wf, hasJacPre);CHKERRQ(ierr); 1472 PetscFunctionReturn(0); 1473 } 1474 1475 /*@C 1476 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. 1477 1478 Not collective 1479 1480 Input Parameters: 1481 + ds - The PetscDS 1482 . f - The test field number 1483 - g - The field number 1484 1485 Output Parameters: 1486 + g0 - integrand for the test and basis function term 1487 . g1 - integrand for the test function and basis function gradient term 1488 . g2 - integrand for the test function gradient and basis function term 1489 - g3 - integrand for the test function gradient and basis function gradient term 1490 1491 Note: We are using a first order FEM model for the weak form: 1492 1493 \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 1494 1495 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1496 1497 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1498 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1499 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1500 $ PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]) 1501 1502 + dim - the spatial dimension 1503 . Nf - the number of fields 1504 . uOff - the offset into u[] and u_t[] for each field 1505 . uOff_x - the offset into u_x[] for each field 1506 . u - each field evaluated at the current point 1507 . u_t - the time derivative of each field evaluated at the current point 1508 . u_x - the gradient of each field evaluated at the current point 1509 . aOff - the offset into a[] and a_t[] for each auxiliary field 1510 . aOff_x - the offset into a_x[] for each auxiliary field 1511 . a - each auxiliary field evaluated at the current point 1512 . a_t - the time derivative of each auxiliary field evaluated at the current point 1513 . a_x - the gradient of auxiliary each field evaluated at the current point 1514 . t - current time 1515 . u_tShift - the multiplier a for dF/dU_t 1516 . x - coordinates of the current point 1517 . numConstants - number of constant parameters 1518 . constants - constant parameters 1519 - g0 - output values at the current point 1520 1521 Level: intermediate 1522 1523 .seealso: PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian() 1524 @*/ 1525 PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, 1526 void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1527 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1528 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1529 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1530 void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1531 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1532 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1533 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1534 void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1535 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1536 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1537 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1538 void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1539 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1540 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1541 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1542 { 1543 PetscPointJac *tmp0, *tmp1, *tmp2, *tmp3; 1544 PetscInt n0, n1, n2, n3; 1545 PetscErrorCode ierr; 1546 1547 PetscFunctionBegin; 1548 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1549 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 1550 if ((g < 0) || (g >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, ds->Nf); 1551 ierr = PetscWeakFormGetJacobianPreconditioner(ds->wf, NULL, 0, f, g, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3);CHKERRQ(ierr); 1552 *g0 = tmp0 ? tmp0[0] : NULL; 1553 *g1 = tmp1 ? tmp1[0] : NULL; 1554 *g2 = tmp2 ? tmp2[0] : NULL; 1555 *g3 = tmp3 ? tmp3[0] : NULL; 1556 PetscFunctionReturn(0); 1557 } 1558 1559 /*@C 1560 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. 1561 1562 Not collective 1563 1564 Input Parameters: 1565 + ds - The PetscDS 1566 . f - The test field number 1567 . g - The field number 1568 . g0 - integrand for the test and basis function term 1569 . g1 - integrand for the test function and basis function gradient term 1570 . g2 - integrand for the test function gradient and basis function term 1571 - g3 - integrand for the test function gradient and basis function gradient term 1572 1573 Note: We are using a first order FEM model for the weak form: 1574 1575 \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 1576 1577 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1578 1579 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1580 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1581 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1582 $ PetscReal t, const PetscReal x[], PetscScalar g0[]) 1583 1584 + dim - the spatial dimension 1585 . Nf - the number of fields 1586 . uOff - the offset into u[] and u_t[] for each field 1587 . uOff_x - the offset into u_x[] for each field 1588 . u - each field evaluated at the current point 1589 . u_t - the time derivative of each field evaluated at the current point 1590 . u_x - the gradient of each field evaluated at the current point 1591 . aOff - the offset into a[] and a_t[] for each auxiliary field 1592 . aOff_x - the offset into a_x[] for each auxiliary field 1593 . a - each auxiliary field evaluated at the current point 1594 . a_t - the time derivative of each auxiliary field evaluated at the current point 1595 . a_x - the gradient of auxiliary each field evaluated at the current point 1596 . t - current time 1597 . u_tShift - the multiplier a for dF/dU_t 1598 . x - coordinates of the current point 1599 . numConstants - number of constant parameters 1600 . constants - constant parameters 1601 - g0 - output values at the current point 1602 1603 Level: intermediate 1604 1605 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobian() 1606 @*/ 1607 PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, 1608 void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1609 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1610 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1611 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1612 void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1613 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1614 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1615 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1616 void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1617 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1618 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1619 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1620 void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1621 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1622 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1623 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1624 { 1625 PetscErrorCode ierr; 1626 1627 PetscFunctionBegin; 1628 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1629 if (g0) PetscValidFunction(g0, 4); 1630 if (g1) PetscValidFunction(g1, 5); 1631 if (g2) PetscValidFunction(g2, 6); 1632 if (g3) PetscValidFunction(g3, 7); 1633 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1634 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 1635 ierr = PetscWeakFormSetIndexJacobianPreconditioner(ds->wf, NULL, 0, f, g, 0, g0, 0, g1, 0, g2, 0, g3);CHKERRQ(ierr); 1636 PetscFunctionReturn(0); 1637 } 1638 1639 /*@C 1640 PetscDSHasDynamicJacobian - Signals that a dynamic Jacobian, dF/du_t, has been set 1641 1642 Not collective 1643 1644 Input Parameter: 1645 . ds - The PetscDS 1646 1647 Output Parameter: 1648 . hasDynJac - flag that pointwise function for dynamic Jacobian has been set 1649 1650 Level: intermediate 1651 1652 .seealso: PetscDSGetDynamicJacobian(), PetscDSSetDynamicJacobian(), PetscDSGetJacobian() 1653 @*/ 1654 PetscErrorCode PetscDSHasDynamicJacobian(PetscDS ds, PetscBool *hasDynJac) 1655 { 1656 PetscErrorCode ierr; 1657 1658 PetscFunctionBegin; 1659 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1660 ierr = PetscWeakFormHasDynamicJacobian(ds->wf, hasDynJac);CHKERRQ(ierr); 1661 PetscFunctionReturn(0); 1662 } 1663 1664 /*@C 1665 PetscDSGetDynamicJacobian - Get the pointwise dynamic Jacobian, dF/du_t, function for given test and basis field 1666 1667 Not collective 1668 1669 Input Parameters: 1670 + ds - The PetscDS 1671 . f - The test field number 1672 - g - The field number 1673 1674 Output Parameters: 1675 + g0 - integrand for the test and basis function term 1676 . g1 - integrand for the test function and basis function gradient term 1677 . g2 - integrand for the test function gradient and basis function term 1678 - g3 - integrand for the test function gradient and basis function gradient term 1679 1680 Note: We are using a first order FEM model for the weak form: 1681 1682 \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 1683 1684 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1685 1686 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1687 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1688 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1689 $ PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]) 1690 1691 + dim - the spatial dimension 1692 . Nf - the number of fields 1693 . uOff - the offset into u[] and u_t[] for each field 1694 . uOff_x - the offset into u_x[] for each field 1695 . u - each field evaluated at the current point 1696 . u_t - the time derivative of each field evaluated at the current point 1697 . u_x - the gradient of each field evaluated at the current point 1698 . aOff - the offset into a[] and a_t[] for each auxiliary field 1699 . aOff_x - the offset into a_x[] for each auxiliary field 1700 . a - each auxiliary field evaluated at the current point 1701 . a_t - the time derivative of each auxiliary field evaluated at the current point 1702 . a_x - the gradient of auxiliary each field evaluated at the current point 1703 . t - current time 1704 . u_tShift - the multiplier a for dF/dU_t 1705 . x - coordinates of the current point 1706 . numConstants - number of constant parameters 1707 . constants - constant parameters 1708 - g0 - output values at the current point 1709 1710 Level: intermediate 1711 1712 .seealso: PetscDSSetJacobian() 1713 @*/ 1714 PetscErrorCode PetscDSGetDynamicJacobian(PetscDS ds, PetscInt f, PetscInt g, 1715 void (**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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1719 void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1720 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1721 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1722 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1723 void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1724 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1725 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1726 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1727 void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1728 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1729 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1730 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1731 { 1732 PetscPointJac *tmp0, *tmp1, *tmp2, *tmp3; 1733 PetscInt n0, n1, n2, n3; 1734 PetscErrorCode ierr; 1735 1736 PetscFunctionBegin; 1737 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1738 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 1739 if ((g < 0) || (g >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, ds->Nf); 1740 ierr = PetscWeakFormGetDynamicJacobian(ds->wf, NULL, 0, f, g, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3);CHKERRQ(ierr); 1741 *g0 = tmp0 ? tmp0[0] : NULL; 1742 *g1 = tmp1 ? tmp1[0] : NULL; 1743 *g2 = tmp2 ? tmp2[0] : NULL; 1744 *g3 = tmp3 ? tmp3[0] : NULL; 1745 PetscFunctionReturn(0); 1746 } 1747 1748 /*@C 1749 PetscDSSetDynamicJacobian - Set the pointwise dynamic Jacobian, dF/du_t, function for given test and basis fields 1750 1751 Not collective 1752 1753 Input Parameters: 1754 + ds - The PetscDS 1755 . f - The test field number 1756 . g - The field number 1757 . g0 - integrand for the test and basis function term 1758 . g1 - integrand for the test function and basis function gradient term 1759 . g2 - integrand for the test function gradient and basis function term 1760 - g3 - integrand for the test function gradient and basis function gradient term 1761 1762 Note: We are using a first order FEM model for the weak form: 1763 1764 \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 1765 1766 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 1767 1768 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1769 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1770 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1771 $ PetscReal t, const PetscReal x[], PetscScalar g0[]) 1772 1773 + dim - the spatial dimension 1774 . Nf - the number of fields 1775 . uOff - the offset into u[] and u_t[] for each field 1776 . uOff_x - the offset into u_x[] for each field 1777 . u - each field evaluated at the current point 1778 . u_t - the time derivative of each field evaluated at the current point 1779 . u_x - the gradient of each field evaluated at the current point 1780 . aOff - the offset into a[] and a_t[] for each auxiliary field 1781 . aOff_x - the offset into a_x[] for each auxiliary field 1782 . a - each auxiliary field evaluated at the current point 1783 . a_t - the time derivative of each auxiliary field evaluated at the current point 1784 . a_x - the gradient of auxiliary each field evaluated at the current point 1785 . t - current time 1786 . u_tShift - the multiplier a for dF/dU_t 1787 . x - coordinates of the current point 1788 . numConstants - number of constant parameters 1789 . constants - constant parameters 1790 - g0 - output values at the current point 1791 1792 Level: intermediate 1793 1794 .seealso: PetscDSGetJacobian() 1795 @*/ 1796 PetscErrorCode PetscDSSetDynamicJacobian(PetscDS ds, PetscInt f, PetscInt g, 1797 void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1798 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1799 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1800 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 1801 void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1802 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1803 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1804 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 1805 void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1806 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1807 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1808 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 1809 void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1810 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1811 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1812 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 1813 { 1814 PetscErrorCode ierr; 1815 1816 PetscFunctionBegin; 1817 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1818 if (g0) PetscValidFunction(g0, 4); 1819 if (g1) PetscValidFunction(g1, 5); 1820 if (g2) PetscValidFunction(g2, 6); 1821 if (g3) PetscValidFunction(g3, 7); 1822 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1823 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 1824 ierr = PetscWeakFormSetIndexDynamicJacobian(ds->wf, NULL, 0, f, g, 0, g0, 0, g1, 0, g2, 0, g3);CHKERRQ(ierr); 1825 PetscFunctionReturn(0); 1826 } 1827 1828 /*@C 1829 PetscDSGetRiemannSolver - Returns the Riemann solver for the given field 1830 1831 Not collective 1832 1833 Input Arguments: 1834 + ds - The PetscDS object 1835 - f - The field number 1836 1837 Output Argument: 1838 . r - Riemann solver 1839 1840 Calling sequence for r: 1841 1842 $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 1843 1844 + dim - The spatial dimension 1845 . Nf - The number of fields 1846 . x - The coordinates at a point on the interface 1847 . n - The normal vector to the interface 1848 . uL - The state vector to the left of the interface 1849 . uR - The state vector to the right of the interface 1850 . flux - output array of flux through the interface 1851 . numConstants - number of constant parameters 1852 . constants - constant parameters 1853 - ctx - optional user context 1854 1855 Level: intermediate 1856 1857 .seealso: PetscDSSetRiemannSolver() 1858 @*/ 1859 PetscErrorCode PetscDSGetRiemannSolver(PetscDS ds, PetscInt f, 1860 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)) 1861 { 1862 PetscRiemannFunc *tmp; 1863 PetscInt n; 1864 PetscErrorCode ierr; 1865 1866 PetscFunctionBegin; 1867 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1868 PetscValidPointer(r, 3); 1869 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 1870 ierr = PetscWeakFormGetRiemannSolver(ds->wf, NULL, 0, f, &n, &tmp);CHKERRQ(ierr); 1871 *r = tmp ? tmp[0] : NULL; 1872 PetscFunctionReturn(0); 1873 } 1874 1875 /*@C 1876 PetscDSSetRiemannSolver - Sets the Riemann solver for the given field 1877 1878 Not collective 1879 1880 Input Arguments: 1881 + ds - The PetscDS object 1882 . f - The field number 1883 - r - Riemann solver 1884 1885 Calling sequence for r: 1886 1887 $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 1888 1889 + dim - The spatial dimension 1890 . Nf - The number of fields 1891 . x - The coordinates at a point on the interface 1892 . n - The normal vector to the interface 1893 . uL - The state vector to the left of the interface 1894 . uR - The state vector to the right of the interface 1895 . flux - output array of flux through the interface 1896 . numConstants - number of constant parameters 1897 . constants - constant parameters 1898 - ctx - optional user context 1899 1900 Level: intermediate 1901 1902 .seealso: PetscDSGetRiemannSolver() 1903 @*/ 1904 PetscErrorCode PetscDSSetRiemannSolver(PetscDS ds, PetscInt f, 1905 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)) 1906 { 1907 PetscErrorCode ierr; 1908 1909 PetscFunctionBegin; 1910 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1911 if (r) PetscValidFunction(r, 3); 1912 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 1913 ierr = PetscWeakFormSetIndexRiemannSolver(ds->wf, NULL, 0, f, 0, r);CHKERRQ(ierr); 1914 PetscFunctionReturn(0); 1915 } 1916 1917 /*@C 1918 PetscDSGetUpdate - Get the pointwise update function for a given field 1919 1920 Not collective 1921 1922 Input Parameters: 1923 + ds - The PetscDS 1924 - f - The field number 1925 1926 Output Parameters: 1927 . update - update function 1928 1929 Note: The calling sequence for the callback update is given by: 1930 1931 $ update(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1932 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1933 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1934 $ PetscReal t, const PetscReal x[], PetscScalar uNew[]) 1935 1936 + dim - the spatial dimension 1937 . Nf - the number of fields 1938 . uOff - the offset into u[] and u_t[] for each field 1939 . uOff_x - the offset into u_x[] for each field 1940 . u - each field evaluated at the current point 1941 . u_t - the time derivative of each field evaluated at the current point 1942 . u_x - the gradient of each field evaluated at the current point 1943 . aOff - the offset into a[] and a_t[] for each auxiliary field 1944 . aOff_x - the offset into a_x[] for each auxiliary field 1945 . a - each auxiliary field evaluated at the current point 1946 . a_t - the time derivative of each auxiliary field evaluated at the current point 1947 . a_x - the gradient of auxiliary each field evaluated at the current point 1948 . t - current time 1949 . x - coordinates of the current point 1950 - uNew - new value for field at the current point 1951 1952 Level: intermediate 1953 1954 .seealso: PetscDSSetUpdate(), PetscDSSetResidual() 1955 @*/ 1956 PetscErrorCode PetscDSGetUpdate(PetscDS ds, PetscInt f, 1957 void (**update)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1958 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1959 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1960 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uNew[])) 1961 { 1962 PetscFunctionBegin; 1963 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1964 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 1965 if (update) {PetscValidPointer(update, 3); *update = ds->update[f];} 1966 PetscFunctionReturn(0); 1967 } 1968 1969 /*@C 1970 PetscDSSetUpdate - Set the pointwise update function for a given field 1971 1972 Not collective 1973 1974 Input Parameters: 1975 + ds - The PetscDS 1976 . f - The field number 1977 - update - update function 1978 1979 Note: The calling sequence for the callback update is given by: 1980 1981 $ update(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1982 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1983 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1984 $ PetscReal t, const PetscReal x[], PetscScalar uNew[]) 1985 1986 + dim - the spatial dimension 1987 . Nf - the number of fields 1988 . uOff - the offset into u[] and u_t[] for each field 1989 . uOff_x - the offset into u_x[] for each field 1990 . u - each field evaluated at the current point 1991 . u_t - the time derivative of each field evaluated at the current point 1992 . u_x - the gradient of each field evaluated at the current point 1993 . aOff - the offset into a[] and a_t[] for each auxiliary field 1994 . aOff_x - the offset into a_x[] for each auxiliary field 1995 . a - each auxiliary field evaluated at the current point 1996 . a_t - the time derivative of each auxiliary field evaluated at the current point 1997 . a_x - the gradient of auxiliary each field evaluated at the current point 1998 . t - current time 1999 . x - coordinates of the current point 2000 - uNew - new field values at the current point 2001 2002 Level: intermediate 2003 2004 .seealso: PetscDSGetResidual() 2005 @*/ 2006 PetscErrorCode PetscDSSetUpdate(PetscDS ds, PetscInt f, 2007 void (*update)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2008 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2009 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2010 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uNew[])) 2011 { 2012 PetscErrorCode ierr; 2013 2014 PetscFunctionBegin; 2015 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2016 if (update) PetscValidFunction(update, 3); 2017 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 2018 ierr = PetscDSEnlarge_Static(ds, f+1);CHKERRQ(ierr); 2019 ds->update[f] = update; 2020 PetscFunctionReturn(0); 2021 } 2022 2023 PetscErrorCode PetscDSGetContext(PetscDS ds, PetscInt f, void **ctx) 2024 { 2025 PetscFunctionBegin; 2026 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2027 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 2028 PetscValidPointer(ctx, 3); 2029 *ctx = ds->ctx[f]; 2030 PetscFunctionReturn(0); 2031 } 2032 2033 PetscErrorCode PetscDSSetContext(PetscDS ds, PetscInt f, void *ctx) 2034 { 2035 PetscErrorCode ierr; 2036 2037 PetscFunctionBegin; 2038 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2039 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 2040 ierr = PetscDSEnlarge_Static(ds, f+1);CHKERRQ(ierr); 2041 ds->ctx[f] = ctx; 2042 PetscFunctionReturn(0); 2043 } 2044 2045 /*@C 2046 PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field 2047 2048 Not collective 2049 2050 Input Parameters: 2051 + ds - The PetscDS 2052 - f - The test field number 2053 2054 Output Parameters: 2055 + f0 - boundary integrand for the test function term 2056 - f1 - boundary integrand for the test function gradient term 2057 2058 Note: We are using a first order FEM model for the weak form: 2059 2060 \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 2061 2062 The calling sequence for the callbacks f0 and f1 is given by: 2063 2064 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2065 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2066 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2067 $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]) 2068 2069 + dim - the spatial dimension 2070 . Nf - the number of fields 2071 . uOff - the offset into u[] and u_t[] for each field 2072 . uOff_x - the offset into u_x[] for each field 2073 . u - each field evaluated at the current point 2074 . u_t - the time derivative of each field evaluated at the current point 2075 . u_x - the gradient of each field evaluated at the current point 2076 . aOff - the offset into a[] and a_t[] for each auxiliary field 2077 . aOff_x - the offset into a_x[] for each auxiliary field 2078 . a - each auxiliary field evaluated at the current point 2079 . a_t - the time derivative of each auxiliary field evaluated at the current point 2080 . a_x - the gradient of auxiliary each field evaluated at the current point 2081 . t - current time 2082 . x - coordinates of the current point 2083 . n - unit normal at the current point 2084 . numConstants - number of constant parameters 2085 . constants - constant parameters 2086 - f0 - output values at the current point 2087 2088 Level: intermediate 2089 2090 .seealso: PetscDSSetBdResidual() 2091 @*/ 2092 PetscErrorCode PetscDSGetBdResidual(PetscDS ds, PetscInt f, 2093 void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2094 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2095 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2096 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), 2097 void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2098 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2099 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2100 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 2101 { 2102 PetscBdPointFunc *tmp0, *tmp1; 2103 PetscInt n0, n1; 2104 PetscErrorCode ierr; 2105 2106 PetscFunctionBegin; 2107 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2108 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 2109 ierr = PetscWeakFormGetBdResidual(ds->wf, NULL, 0, f, &n0, &tmp0, &n1, &tmp1);CHKERRQ(ierr); 2110 *f0 = tmp0 ? tmp0[0] : NULL; 2111 *f1 = tmp1 ? tmp1[0] : NULL; 2112 PetscFunctionReturn(0); 2113 } 2114 2115 /*@C 2116 PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field 2117 2118 Not collective 2119 2120 Input Parameters: 2121 + ds - The PetscDS 2122 . f - The test field number 2123 . f0 - boundary integrand for the test function term 2124 - f1 - boundary integrand for the test function gradient term 2125 2126 Note: We are using a first order FEM model for the weak form: 2127 2128 \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 2129 2130 The calling sequence for the callbacks f0 and f1 is given by: 2131 2132 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2133 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2134 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2135 $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]) 2136 2137 + dim - the spatial dimension 2138 . Nf - the number of fields 2139 . uOff - the offset into u[] and u_t[] for each field 2140 . uOff_x - the offset into u_x[] for each field 2141 . u - each field evaluated at the current point 2142 . u_t - the time derivative of each field evaluated at the current point 2143 . u_x - the gradient of each field evaluated at the current point 2144 . aOff - the offset into a[] and a_t[] for each auxiliary field 2145 . aOff_x - the offset into a_x[] for each auxiliary field 2146 . a - each auxiliary field evaluated at the current point 2147 . a_t - the time derivative of each auxiliary field evaluated at the current point 2148 . a_x - the gradient of auxiliary each field evaluated at the current point 2149 . t - current time 2150 . x - coordinates of the current point 2151 . n - unit normal at the current point 2152 . numConstants - number of constant parameters 2153 . constants - constant parameters 2154 - f0 - output values at the current point 2155 2156 Level: intermediate 2157 2158 .seealso: PetscDSGetBdResidual() 2159 @*/ 2160 PetscErrorCode PetscDSSetBdResidual(PetscDS ds, PetscInt f, 2161 void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2162 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2163 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2164 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), 2165 void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2166 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2167 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2168 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])) 2169 { 2170 PetscErrorCode ierr; 2171 2172 PetscFunctionBegin; 2173 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2174 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 2175 ierr = PetscWeakFormSetIndexBdResidual(ds->wf, NULL, 0, f, 0, f0, 0, f1);CHKERRQ(ierr); 2176 PetscFunctionReturn(0); 2177 } 2178 2179 /*@ 2180 PetscDSHasBdJacobian - Signals that boundary Jacobian functions have been set 2181 2182 Not collective 2183 2184 Input Parameter: 2185 . ds - The PetscDS 2186 2187 Output Parameter: 2188 . hasBdJac - flag that pointwise function for the boundary Jacobian has been set 2189 2190 Level: intermediate 2191 2192 .seealso: PetscDSHasJacobian(), PetscDSSetBdJacobian(), PetscDSGetBdJacobian() 2193 @*/ 2194 PetscErrorCode PetscDSHasBdJacobian(PetscDS ds, PetscBool *hasBdJac) 2195 { 2196 PetscErrorCode ierr; 2197 2198 PetscFunctionBegin; 2199 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2200 PetscValidBoolPointer(hasBdJac, 2); 2201 ierr = PetscWeakFormHasBdJacobian(ds->wf, hasBdJac);CHKERRQ(ierr); 2202 PetscFunctionReturn(0); 2203 } 2204 2205 /*@C 2206 PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field 2207 2208 Not collective 2209 2210 Input Parameters: 2211 + ds - The PetscDS 2212 . f - The test field number 2213 - g - The field number 2214 2215 Output Parameters: 2216 + g0 - integrand for the test and basis function term 2217 . g1 - integrand for the test function and basis function gradient term 2218 . g2 - integrand for the test function gradient and basis function term 2219 - g3 - integrand for the test function gradient and basis function gradient term 2220 2221 Note: We are using a first order FEM model for the weak form: 2222 2223 \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 2224 2225 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 2226 2227 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2228 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2229 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2230 $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[]) 2231 2232 + dim - the spatial dimension 2233 . Nf - the number of fields 2234 . uOff - the offset into u[] and u_t[] for each field 2235 . uOff_x - the offset into u_x[] for each field 2236 . u - each field evaluated at the current point 2237 . u_t - the time derivative of each field evaluated at the current point 2238 . u_x - the gradient of each field evaluated at the current point 2239 . aOff - the offset into a[] and a_t[] for each auxiliary field 2240 . aOff_x - the offset into a_x[] for each auxiliary field 2241 . a - each auxiliary field evaluated at the current point 2242 . a_t - the time derivative of each auxiliary field evaluated at the current point 2243 . a_x - the gradient of auxiliary each field evaluated at the current point 2244 . t - current time 2245 . u_tShift - the multiplier a for dF/dU_t 2246 . x - coordinates of the current point 2247 . n - normal at the current point 2248 . numConstants - number of constant parameters 2249 . constants - constant parameters 2250 - g0 - output values at the current point 2251 2252 Level: intermediate 2253 2254 .seealso: PetscDSSetBdJacobian() 2255 @*/ 2256 PetscErrorCode PetscDSGetBdJacobian(PetscDS ds, PetscInt f, PetscInt g, 2257 void (**g0)(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 g0[]), 2261 void (**g1)(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 g1[]), 2265 void (**g2)(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 g2[]), 2269 void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2270 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2271 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2272 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 2273 { 2274 PetscBdPointJac *tmp0, *tmp1, *tmp2, *tmp3; 2275 PetscInt n0, n1, n2, n3; 2276 PetscErrorCode ierr; 2277 2278 PetscFunctionBegin; 2279 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2280 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 2281 if ((g < 0) || (g >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, ds->Nf); 2282 ierr = PetscWeakFormGetBdJacobian(ds->wf, NULL, 0, f, g, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3);CHKERRQ(ierr); 2283 *g0 = tmp0 ? tmp0[0] : NULL; 2284 *g1 = tmp1 ? tmp1[0] : NULL; 2285 *g2 = tmp2 ? tmp2[0] : NULL; 2286 *g3 = tmp3 ? tmp3[0] : NULL; 2287 PetscFunctionReturn(0); 2288 } 2289 2290 /*@C 2291 PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field 2292 2293 Not collective 2294 2295 Input Parameters: 2296 + ds - The PetscDS 2297 . f - The test field number 2298 . g - The field number 2299 . g0 - integrand for the test and basis function term 2300 . g1 - integrand for the test function and basis function gradient term 2301 . g2 - integrand for the test function gradient and basis function term 2302 - g3 - integrand for the test function gradient and basis function gradient term 2303 2304 Note: We are using a first order FEM model for the weak form: 2305 2306 \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 2307 2308 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 2309 2310 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2311 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2312 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2313 $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[]) 2314 2315 + dim - the spatial dimension 2316 . Nf - the number of fields 2317 . uOff - the offset into u[] and u_t[] for each field 2318 . uOff_x - the offset into u_x[] for each field 2319 . u - each field evaluated at the current point 2320 . u_t - the time derivative of each field evaluated at the current point 2321 . u_x - the gradient of each field evaluated at the current point 2322 . aOff - the offset into a[] and a_t[] for each auxiliary field 2323 . aOff_x - the offset into a_x[] for each auxiliary field 2324 . a - each auxiliary field evaluated at the current point 2325 . a_t - the time derivative of each auxiliary field evaluated at the current point 2326 . a_x - the gradient of auxiliary each field evaluated at the current point 2327 . t - current time 2328 . u_tShift - the multiplier a for dF/dU_t 2329 . x - coordinates of the current point 2330 . n - normal at the current point 2331 . numConstants - number of constant parameters 2332 . constants - constant parameters 2333 - g0 - output values at the current point 2334 2335 Level: intermediate 2336 2337 .seealso: PetscDSGetBdJacobian() 2338 @*/ 2339 PetscErrorCode PetscDSSetBdJacobian(PetscDS ds, PetscInt f, PetscInt g, 2340 void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2341 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2342 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2343 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 2344 void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2345 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2346 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2347 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 2348 void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2349 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2350 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2351 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 2352 void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2353 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2354 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2355 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 2356 { 2357 PetscErrorCode ierr; 2358 2359 PetscFunctionBegin; 2360 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2361 if (g0) PetscValidFunction(g0, 4); 2362 if (g1) PetscValidFunction(g1, 5); 2363 if (g2) PetscValidFunction(g2, 6); 2364 if (g3) PetscValidFunction(g3, 7); 2365 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 2366 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 2367 ierr = PetscWeakFormSetIndexBdJacobian(ds->wf, NULL, 0, f, g, 0, g0, 0, g1, 0, g2, 0, g3);CHKERRQ(ierr); 2368 PetscFunctionReturn(0); 2369 } 2370 2371 /*@ 2372 PetscDSHasBdJacobianPreconditioner - Signals that boundary Jacobian preconditioner functions have been set 2373 2374 Not collective 2375 2376 Input Parameter: 2377 . ds - The PetscDS 2378 2379 Output Parameter: 2380 . hasBdJac - flag that pointwise function for the boundary Jacobian preconditioner has been set 2381 2382 Level: intermediate 2383 2384 .seealso: PetscDSHasJacobian(), PetscDSSetBdJacobian(), PetscDSGetBdJacobian() 2385 @*/ 2386 PetscErrorCode PetscDSHasBdJacobianPreconditioner(PetscDS ds, PetscBool *hasBdJacPre) 2387 { 2388 PetscErrorCode ierr; 2389 2390 PetscFunctionBegin; 2391 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2392 PetscValidBoolPointer(hasBdJacPre, 2); 2393 ierr = PetscWeakFormHasBdJacobianPreconditioner(ds->wf, hasBdJacPre);CHKERRQ(ierr); 2394 PetscFunctionReturn(0); 2395 } 2396 2397 /*@C 2398 PetscDSGetBdJacobianPreconditioner - Get the pointwise boundary Jacobian preconditioner function for given test and basis field 2399 2400 Not collective 2401 2402 Input Parameters: 2403 + ds - The PetscDS 2404 . f - The test field number 2405 - g - The field number 2406 2407 Output Parameters: 2408 + g0 - integrand for the test and basis function term 2409 . g1 - integrand for the test function and basis function gradient term 2410 . g2 - integrand for the test function gradient and basis function term 2411 - g3 - integrand for the test function gradient and basis function gradient term 2412 2413 Note: We are using a first order FEM model for the weak form: 2414 2415 \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 2416 2417 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 2418 2419 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2420 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2421 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2422 $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 2423 2424 + dim - the spatial dimension 2425 . Nf - the number of fields 2426 . NfAux - the number of auxiliary fields 2427 . uOff - the offset into u[] and u_t[] for each field 2428 . uOff_x - the offset into u_x[] for each field 2429 . u - each field evaluated at the current point 2430 . u_t - the time derivative of each field evaluated at the current point 2431 . u_x - the gradient of each field evaluated at the current point 2432 . aOff - the offset into a[] and a_t[] for each auxiliary field 2433 . aOff_x - the offset into a_x[] for each auxiliary field 2434 . a - each auxiliary field evaluated at the current point 2435 . a_t - the time derivative of each auxiliary field evaluated at the current point 2436 . a_x - the gradient of auxiliary each field evaluated at the current point 2437 . t - current time 2438 . u_tShift - the multiplier a for dF/dU_t 2439 . x - coordinates of the current point 2440 . n - normal at the current point 2441 . numConstants - number of constant parameters 2442 . constants - constant parameters 2443 - g0 - output values at the current point 2444 2445 This is not yet available in Fortran. 2446 2447 Level: intermediate 2448 2449 .seealso: PetscDSSetBdJacobianPreconditioner() 2450 @*/ 2451 PetscErrorCode PetscDSGetBdJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, 2452 void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2453 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2454 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2455 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 2456 void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2457 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2458 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2459 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 2460 void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2461 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2462 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2463 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 2464 void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2465 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2466 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2467 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 2468 { 2469 PetscBdPointJac *tmp0, *tmp1, *tmp2, *tmp3; 2470 PetscInt n0, n1, n2, n3; 2471 PetscErrorCode ierr; 2472 2473 PetscFunctionBegin; 2474 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2475 if ((f < 0) || (f >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, ds->Nf); 2476 if ((g < 0) || (g >= ds->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, ds->Nf); 2477 ierr = PetscWeakFormGetBdJacobianPreconditioner(ds->wf, NULL, 0, f, g, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3);CHKERRQ(ierr); 2478 *g0 = tmp0 ? tmp0[0] : NULL; 2479 *g1 = tmp1 ? tmp1[0] : NULL; 2480 *g2 = tmp2 ? tmp2[0] : NULL; 2481 *g3 = tmp3 ? tmp3[0] : NULL; 2482 PetscFunctionReturn(0); 2483 } 2484 2485 /*@C 2486 PetscDSSetBdJacobianPreconditioner - Set the pointwise boundary Jacobian preconditioner function for given test and basis field 2487 2488 Not collective 2489 2490 Input Parameters: 2491 + ds - The PetscDS 2492 . f - The test field number 2493 . g - The field number 2494 . g0 - integrand for the test and basis function term 2495 . g1 - integrand for the test function and basis function gradient term 2496 . g2 - integrand for the test function gradient and basis function term 2497 - g3 - integrand for the test function gradient and basis function gradient term 2498 2499 Note: We are using a first order FEM model for the weak form: 2500 2501 \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 2502 2503 The calling sequence for the callbacks g0, g1, g2 and g3 is given by: 2504 2505 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2506 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2507 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2508 $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 2509 2510 + dim - the spatial dimension 2511 . Nf - the number of fields 2512 . NfAux - the number of auxiliary fields 2513 . uOff - the offset into u[] and u_t[] for each field 2514 . uOff_x - the offset into u_x[] for each field 2515 . u - each field evaluated at the current point 2516 . u_t - the time derivative of each field evaluated at the current point 2517 . u_x - the gradient of each field evaluated at the current point 2518 . aOff - the offset into a[] and a_t[] for each auxiliary field 2519 . aOff_x - the offset into a_x[] for each auxiliary field 2520 . a - each auxiliary field evaluated at the current point 2521 . a_t - the time derivative of each auxiliary field evaluated at the current point 2522 . a_x - the gradient of auxiliary each field evaluated at the current point 2523 . t - current time 2524 . u_tShift - the multiplier a for dF/dU_t 2525 . x - coordinates of the current point 2526 . n - normal at the current point 2527 . numConstants - number of constant parameters 2528 . constants - constant parameters 2529 - g0 - output values at the current point 2530 2531 This is not yet available in Fortran. 2532 2533 Level: intermediate 2534 2535 .seealso: PetscDSGetBdJacobianPreconditioner() 2536 @*/ 2537 PetscErrorCode PetscDSSetBdJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, 2538 void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2539 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2540 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2541 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), 2542 void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2543 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2544 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2545 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]), 2546 void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2547 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2548 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2549 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]), 2550 void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2551 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2552 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2553 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])) 2554 { 2555 PetscErrorCode ierr; 2556 2557 PetscFunctionBegin; 2558 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 2559 if (g0) PetscValidFunction(g0, 4); 2560 if (g1) PetscValidFunction(g1, 5); 2561 if (g2) PetscValidFunction(g2, 6); 2562 if (g3) PetscValidFunction(g3, 7); 2563 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 2564 if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g); 2565 ierr = PetscWeakFormSetIndexBdJacobianPreconditioner(ds->wf, NULL, 0, f, g, 0, g0, 0, g1, 0, g2, 0, g3);CHKERRQ(ierr); 2566 PetscFunctionReturn(0); 2567 } 2568 2569 /*@C 2570 PetscDSGetExactSolution - Get the pointwise exact solution function for a given test field 2571 2572 Not collective 2573 2574 Input Parameters: 2575 + prob - The PetscDS 2576 - f - The test field number 2577 2578 Output Parameter: 2579 + exactSol - exact solution for the test field 2580 - exactCtx - exact solution context 2581 2582 Note: The calling sequence for the solution functions is given by: 2583 2584 $ sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) 2585 2586 + dim - the spatial dimension 2587 . t - current time 2588 . x - coordinates of the current point 2589 . Nc - the number of field components 2590 . u - the solution field evaluated at the current point 2591 - ctx - a user context 2592 2593 Level: intermediate 2594 2595 .seealso: PetscDSSetExactSolution(), PetscDSGetExactSolutionTimeDerivative() 2596 @*/ 2597 PetscErrorCode PetscDSGetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (**sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void **ctx) 2598 { 2599 PetscFunctionBegin; 2600 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2601 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); 2602 if (sol) {PetscValidPointer(sol, 3); *sol = prob->exactSol[f];} 2603 if (ctx) {PetscValidPointer(ctx, 4); *ctx = prob->exactCtx[f];} 2604 PetscFunctionReturn(0); 2605 } 2606 2607 /*@C 2608 PetscDSSetExactSolution - Set the pointwise exact solution function for a given test field 2609 2610 Not collective 2611 2612 Input Parameters: 2613 + prob - The PetscDS 2614 . f - The test field number 2615 . sol - solution function for the test fields 2616 - ctx - solution context or NULL 2617 2618 Note: The calling sequence for solution functions is given by: 2619 2620 $ sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) 2621 2622 + dim - the spatial dimension 2623 . t - current time 2624 . x - coordinates of the current point 2625 . Nc - the number of field components 2626 . u - the solution field evaluated at the current point 2627 - ctx - a user context 2628 2629 Level: intermediate 2630 2631 .seealso: PetscDSGetExactSolution() 2632 @*/ 2633 PetscErrorCode PetscDSSetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (*sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void *ctx) 2634 { 2635 PetscErrorCode ierr; 2636 2637 PetscFunctionBegin; 2638 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2639 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 2640 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 2641 if (sol) {PetscValidFunction(sol, 3); prob->exactSol[f] = sol;} 2642 if (ctx) {PetscValidFunction(ctx, 4); prob->exactCtx[f] = ctx;} 2643 PetscFunctionReturn(0); 2644 } 2645 2646 /*@C 2647 PetscDSGetExactSolutionTimeDerivative - Get the pointwise time derivative of the exact solution function for a given test field 2648 2649 Not collective 2650 2651 Input Parameters: 2652 + prob - The PetscDS 2653 - f - The test field number 2654 2655 Output Parameter: 2656 + exactSol - time derivative of the exact solution for the test field 2657 - exactCtx - time derivative of the exact solution context 2658 2659 Note: The calling sequence for the solution functions is given by: 2660 2661 $ sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) 2662 2663 + dim - the spatial dimension 2664 . t - current time 2665 . x - coordinates of the current point 2666 . Nc - the number of field components 2667 . u - the solution field evaluated at the current point 2668 - ctx - a user context 2669 2670 Level: intermediate 2671 2672 .seealso: PetscDSSetExactSolutionTimeDerivative(), PetscDSGetExactSolution() 2673 @*/ 2674 PetscErrorCode PetscDSGetExactSolutionTimeDerivative(PetscDS prob, PetscInt f, PetscErrorCode (**sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void **ctx) 2675 { 2676 PetscFunctionBegin; 2677 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2678 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); 2679 if (sol) {PetscValidPointer(sol, 3); *sol = prob->exactSol_t[f];} 2680 if (ctx) {PetscValidPointer(ctx, 4); *ctx = prob->exactCtx_t[f];} 2681 PetscFunctionReturn(0); 2682 } 2683 2684 /*@C 2685 PetscDSSetExactSolutionTimeDerivative - Set the pointwise time derivative of the exact solution function for a given test field 2686 2687 Not collective 2688 2689 Input Parameters: 2690 + prob - The PetscDS 2691 . f - The test field number 2692 . sol - time derivative of the solution function for the test fields 2693 - ctx - time derivative of the solution context or NULL 2694 2695 Note: The calling sequence for solution functions is given by: 2696 2697 $ sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) 2698 2699 + dim - the spatial dimension 2700 . t - current time 2701 . x - coordinates of the current point 2702 . Nc - the number of field components 2703 . u - the solution field evaluated at the current point 2704 - ctx - a user context 2705 2706 Level: intermediate 2707 2708 .seealso: PetscDSGetExactSolutionTimeDerivative(), PetscDSSetExactSolution() 2709 @*/ 2710 PetscErrorCode PetscDSSetExactSolutionTimeDerivative(PetscDS prob, PetscInt f, PetscErrorCode (*sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void *ctx) 2711 { 2712 PetscErrorCode ierr; 2713 2714 PetscFunctionBegin; 2715 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2716 if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f); 2717 ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr); 2718 if (sol) {PetscValidFunction(sol, 3); prob->exactSol_t[f] = sol;} 2719 if (ctx) {PetscValidFunction(ctx, 4); prob->exactCtx_t[f] = ctx;} 2720 PetscFunctionReturn(0); 2721 } 2722 2723 /*@C 2724 PetscDSGetConstants - Returns the array of constants passed to point functions 2725 2726 Not collective 2727 2728 Input Parameter: 2729 . prob - The PetscDS object 2730 2731 Output Parameters: 2732 + numConstants - The number of constants 2733 - constants - The array of constants, NULL if there are none 2734 2735 Level: intermediate 2736 2737 .seealso: PetscDSSetConstants(), PetscDSCreate() 2738 @*/ 2739 PetscErrorCode PetscDSGetConstants(PetscDS prob, PetscInt *numConstants, const PetscScalar *constants[]) 2740 { 2741 PetscFunctionBegin; 2742 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2743 if (numConstants) {PetscValidPointer(numConstants, 2); *numConstants = prob->numConstants;} 2744 if (constants) {PetscValidPointer(constants, 3); *constants = prob->constants;} 2745 PetscFunctionReturn(0); 2746 } 2747 2748 /*@C 2749 PetscDSSetConstants - Set the array of constants passed to point functions 2750 2751 Not collective 2752 2753 Input Parameters: 2754 + prob - The PetscDS object 2755 . numConstants - The number of constants 2756 - constants - The array of constants, NULL if there are none 2757 2758 Level: intermediate 2759 2760 .seealso: PetscDSGetConstants(), PetscDSCreate() 2761 @*/ 2762 PetscErrorCode PetscDSSetConstants(PetscDS prob, PetscInt numConstants, PetscScalar constants[]) 2763 { 2764 PetscErrorCode ierr; 2765 2766 PetscFunctionBegin; 2767 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2768 if (numConstants != prob->numConstants) { 2769 ierr = PetscFree(prob->constants);CHKERRQ(ierr); 2770 prob->numConstants = numConstants; 2771 if (prob->numConstants) { 2772 ierr = PetscMalloc1(prob->numConstants, &prob->constants);CHKERRQ(ierr); 2773 } else { 2774 prob->constants = NULL; 2775 } 2776 } 2777 if (prob->numConstants) { 2778 PetscValidPointer(constants, 3); 2779 ierr = PetscArraycpy(prob->constants, constants, prob->numConstants);CHKERRQ(ierr); 2780 } 2781 PetscFunctionReturn(0); 2782 } 2783 2784 /*@ 2785 PetscDSGetFieldIndex - Returns the index of the given field 2786 2787 Not collective 2788 2789 Input Parameters: 2790 + prob - The PetscDS object 2791 - disc - The discretization object 2792 2793 Output Parameter: 2794 . f - The field number 2795 2796 Level: beginner 2797 2798 .seealso: PetscGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate() 2799 @*/ 2800 PetscErrorCode PetscDSGetFieldIndex(PetscDS prob, PetscObject disc, PetscInt *f) 2801 { 2802 PetscInt g; 2803 2804 PetscFunctionBegin; 2805 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2806 PetscValidPointer(f, 3); 2807 *f = -1; 2808 for (g = 0; g < prob->Nf; ++g) {if (disc == prob->disc[g]) break;} 2809 if (g == prob->Nf) SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Field not found in PetscDS."); 2810 *f = g; 2811 PetscFunctionReturn(0); 2812 } 2813 2814 /*@ 2815 PetscDSGetFieldSize - Returns the size of the given field in the full space basis 2816 2817 Not collective 2818 2819 Input Parameters: 2820 + prob - The PetscDS object 2821 - f - The field number 2822 2823 Output Parameter: 2824 . size - The size 2825 2826 Level: beginner 2827 2828 .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate() 2829 @*/ 2830 PetscErrorCode PetscDSGetFieldSize(PetscDS prob, PetscInt f, PetscInt *size) 2831 { 2832 PetscErrorCode ierr; 2833 2834 PetscFunctionBegin; 2835 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2836 PetscValidPointer(size, 3); 2837 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); 2838 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2839 *size = prob->Nb[f]; 2840 PetscFunctionReturn(0); 2841 } 2842 2843 /*@ 2844 PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis 2845 2846 Not collective 2847 2848 Input Parameters: 2849 + prob - The PetscDS object 2850 - f - The field number 2851 2852 Output Parameter: 2853 . off - The offset 2854 2855 Level: beginner 2856 2857 .seealso: PetscDSGetFieldSize(), PetscDSGetNumFields(), PetscDSCreate() 2858 @*/ 2859 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off) 2860 { 2861 PetscInt size, g; 2862 PetscErrorCode ierr; 2863 2864 PetscFunctionBegin; 2865 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2866 PetscValidPointer(off, 3); 2867 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); 2868 *off = 0; 2869 for (g = 0; g < f; ++g) { 2870 ierr = PetscDSGetFieldSize(prob, g, &size);CHKERRQ(ierr); 2871 *off += size; 2872 } 2873 PetscFunctionReturn(0); 2874 } 2875 2876 /*@ 2877 PetscDSGetDimensions - Returns the size of the approximation space for each field on an evaluation point 2878 2879 Not collective 2880 2881 Input Parameter: 2882 . prob - The PetscDS object 2883 2884 Output Parameter: 2885 . dimensions - The number of dimensions 2886 2887 Level: beginner 2888 2889 .seealso: PetscDSGetComponentOffsets(), PetscDSGetNumFields(), PetscDSCreate() 2890 @*/ 2891 PetscErrorCode PetscDSGetDimensions(PetscDS prob, PetscInt *dimensions[]) 2892 { 2893 PetscErrorCode ierr; 2894 2895 PetscFunctionBegin; 2896 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2897 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2898 PetscValidPointer(dimensions, 2); 2899 *dimensions = prob->Nb; 2900 PetscFunctionReturn(0); 2901 } 2902 2903 /*@ 2904 PetscDSGetComponents - Returns the number of components for each field on an evaluation point 2905 2906 Not collective 2907 2908 Input Parameter: 2909 . prob - The PetscDS object 2910 2911 Output Parameter: 2912 . components - The number of components 2913 2914 Level: beginner 2915 2916 .seealso: PetscDSGetComponentOffsets(), PetscDSGetNumFields(), PetscDSCreate() 2917 @*/ 2918 PetscErrorCode PetscDSGetComponents(PetscDS prob, PetscInt *components[]) 2919 { 2920 PetscErrorCode ierr; 2921 2922 PetscFunctionBegin; 2923 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2924 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 2925 PetscValidPointer(components, 2); 2926 *components = prob->Nc; 2927 PetscFunctionReturn(0); 2928 } 2929 2930 /*@ 2931 PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point 2932 2933 Not collective 2934 2935 Input Parameters: 2936 + prob - The PetscDS object 2937 - f - The field number 2938 2939 Output Parameter: 2940 . off - The offset 2941 2942 Level: beginner 2943 2944 .seealso: PetscDSGetNumFields(), PetscDSCreate() 2945 @*/ 2946 PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off) 2947 { 2948 PetscFunctionBegin; 2949 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2950 PetscValidPointer(off, 3); 2951 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); 2952 *off = prob->off[f]; 2953 PetscFunctionReturn(0); 2954 } 2955 2956 /*@ 2957 PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point 2958 2959 Not collective 2960 2961 Input Parameter: 2962 . prob - The PetscDS object 2963 2964 Output Parameter: 2965 . offsets - The offsets 2966 2967 Level: beginner 2968 2969 .seealso: PetscDSGetNumFields(), PetscDSCreate() 2970 @*/ 2971 PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[]) 2972 { 2973 PetscFunctionBegin; 2974 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2975 PetscValidPointer(offsets, 2); 2976 *offsets = prob->off; 2977 PetscFunctionReturn(0); 2978 } 2979 2980 /*@ 2981 PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point 2982 2983 Not collective 2984 2985 Input Parameter: 2986 . prob - The PetscDS object 2987 2988 Output Parameter: 2989 . offsets - The offsets 2990 2991 Level: beginner 2992 2993 .seealso: PetscDSGetNumFields(), PetscDSCreate() 2994 @*/ 2995 PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[]) 2996 { 2997 PetscFunctionBegin; 2998 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 2999 PetscValidPointer(offsets, 2); 3000 *offsets = prob->offDer; 3001 PetscFunctionReturn(0); 3002 } 3003 3004 /*@C 3005 PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization 3006 3007 Not collective 3008 3009 Input Parameter: 3010 . prob - The PetscDS object 3011 3012 Output Parameter: 3013 . T - The basis function and derivatives tabulation at quadrature points for each field 3014 3015 Level: intermediate 3016 3017 .seealso: PetscDSCreate() 3018 @*/ 3019 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscTabulation *T[]) 3020 { 3021 PetscErrorCode ierr; 3022 3023 PetscFunctionBegin; 3024 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3025 PetscValidPointer(T, 2); 3026 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 3027 *T = prob->T; 3028 PetscFunctionReturn(0); 3029 } 3030 3031 /*@C 3032 PetscDSGetFaceTabulation - Return the basis tabulation at quadrature points on the faces 3033 3034 Not collective 3035 3036 Input Parameter: 3037 . prob - The PetscDS object 3038 3039 Output Parameter: 3040 . Tf - The basis function and derviative tabulation on each local face at quadrature points for each and field 3041 3042 Level: intermediate 3043 3044 .seealso: PetscDSGetTabulation(), PetscDSCreate() 3045 @*/ 3046 PetscErrorCode PetscDSGetFaceTabulation(PetscDS prob, PetscTabulation *Tf[]) 3047 { 3048 PetscErrorCode ierr; 3049 3050 PetscFunctionBegin; 3051 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3052 PetscValidPointer(Tf, 2); 3053 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 3054 *Tf = prob->Tf; 3055 PetscFunctionReturn(0); 3056 } 3057 3058 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x) 3059 { 3060 PetscErrorCode ierr; 3061 3062 PetscFunctionBegin; 3063 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3064 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 3065 if (u) {PetscValidPointer(u, 2); *u = prob->u;} 3066 if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;} 3067 if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;} 3068 PetscFunctionReturn(0); 3069 } 3070 3071 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3) 3072 { 3073 PetscErrorCode ierr; 3074 3075 PetscFunctionBegin; 3076 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3077 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 3078 if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;} 3079 if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;} 3080 if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;} 3081 if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;} 3082 if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;} 3083 if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;} 3084 PetscFunctionReturn(0); 3085 } 3086 3087 PetscErrorCode PetscDSGetWorkspace(PetscDS prob, PetscReal **x, PetscScalar **basisReal, PetscScalar **basisDerReal, PetscScalar **testReal, PetscScalar **testDerReal) 3088 { 3089 PetscErrorCode ierr; 3090 3091 PetscFunctionBegin; 3092 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3093 ierr = PetscDSSetUp(prob);CHKERRQ(ierr); 3094 if (x) {PetscValidPointer(x, 2); *x = prob->x;} 3095 if (basisReal) {PetscValidPointer(basisReal, 3); *basisReal = prob->basisReal;} 3096 if (basisDerReal) {PetscValidPointer(basisDerReal, 4); *basisDerReal = prob->basisDerReal;} 3097 if (testReal) {PetscValidPointer(testReal, 5); *testReal = prob->testReal;} 3098 if (testDerReal) {PetscValidPointer(testDerReal, 6); *testDerReal = prob->testDerReal;} 3099 PetscFunctionReturn(0); 3100 } 3101 3102 /*@C 3103 PetscDSAddBoundary - Add a boundary condition to the model. The pointwise functions are used to provide boundary values for essential boundary conditions. In FEM, they are acting upon by dual basis functionals to generate FEM coefficients which are fixed. Natural boundary conditions signal to PETSc that boundary integrals should be performaed, using the kernels from PetscDSSetBdResidual(). 3104 3105 Collective on ds 3106 3107 Input Parameters: 3108 + ds - The PetscDS object 3109 . type - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann) 3110 . name - The BC name 3111 . label - The label defining constrained points 3112 . Nv - The number of DMLabel values for constrained points 3113 . values - An array of label values for constrained points 3114 . field - The field to constrain 3115 . Nc - The number of constrained field components (0 will constrain all fields) 3116 . comps - An array of constrained component numbers 3117 . bcFunc - A pointwise function giving boundary values 3118 . bcFunc_t - A pointwise function giving the time derviative of the boundary values, or NULL 3119 - ctx - An optional user context for bcFunc 3120 3121 Output Parameters: 3122 - bd - The boundary number 3123 3124 Options Database Keys: 3125 + -bc_<boundary name> <num> - Overrides the boundary ids 3126 - -bc_<boundary name>_comp <num> - Overrides the boundary components 3127 3128 Note: 3129 Both bcFunc abd bcFunc_t will depend on the boundary condition type. If the type if DM_BC_ESSENTIAL, Then the calling sequence is: 3130 3131 $ bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[]) 3132 3133 If the type is DM_BC_ESSENTIAL_FIELD or other _FIELD value, then the calling sequence is: 3134 3135 $ bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, 3136 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 3137 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 3138 $ PetscReal time, const PetscReal x[], PetscScalar bcval[]) 3139 3140 + dim - the spatial dimension 3141 . Nf - the number of fields 3142 . uOff - the offset into u[] and u_t[] for each field 3143 . uOff_x - the offset into u_x[] for each field 3144 . u - each field evaluated at the current point 3145 . u_t - the time derivative of each field evaluated at the current point 3146 . u_x - the gradient of each field evaluated at the current point 3147 . aOff - the offset into a[] and a_t[] for each auxiliary field 3148 . aOff_x - the offset into a_x[] for each auxiliary field 3149 . a - each auxiliary field evaluated at the current point 3150 . a_t - the time derivative of each auxiliary field evaluated at the current point 3151 . a_x - the gradient of auxiliary each field evaluated at the current point 3152 . t - current time 3153 . x - coordinates of the current point 3154 . numConstants - number of constant parameters 3155 . constants - constant parameters 3156 - bcval - output values at the current point 3157 3158 Level: developer 3159 3160 .seealso: PetscDSAddBoundaryByName(), PetscDSGetBoundary(), PetscDSSetResidual(), PetscDSSetBdResidual() 3161 @*/ 3162 PetscErrorCode PetscDSAddBoundary(PetscDS ds, DMBoundaryConditionType type, const char name[], DMLabel label, PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx, PetscInt *bd) 3163 { 3164 DSBoundary head = ds->boundary, b; 3165 PetscInt n = 0; 3166 const char *lname; 3167 PetscErrorCode ierr; 3168 3169 PetscFunctionBegin; 3170 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3171 PetscValidLogicalCollectiveEnum(ds, type, 2); 3172 PetscValidCharPointer(name, 3); 3173 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4); 3174 PetscValidLogicalCollectiveInt(ds, Nv, 5); 3175 PetscValidLogicalCollectiveInt(ds, field, 7); 3176 PetscValidLogicalCollectiveInt(ds, Nc, 8); 3177 ierr = PetscNew(&b);CHKERRQ(ierr); 3178 ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr); 3179 ierr = PetscWeakFormCreate(PETSC_COMM_SELF, &b->wf);CHKERRQ(ierr); 3180 ierr = PetscWeakFormSetNumFields(b->wf, ds->Nf);CHKERRQ(ierr); 3181 ierr = PetscMalloc1(Nv, &b->values);CHKERRQ(ierr); 3182 if (Nv) {ierr = PetscArraycpy(b->values, values, Nv);CHKERRQ(ierr);} 3183 ierr = PetscMalloc1(Nc, &b->comps);CHKERRQ(ierr); 3184 if (Nc) {ierr = PetscArraycpy(b->comps, comps, Nc);CHKERRQ(ierr);} 3185 ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr); 3186 ierr = PetscStrallocpy(lname, (char **) &b->lname);CHKERRQ(ierr); 3187 b->type = type; 3188 b->label = label; 3189 b->Nv = Nv; 3190 b->field = field; 3191 b->Nc = Nc; 3192 b->func = bcFunc; 3193 b->func_t = bcFunc_t; 3194 b->ctx = ctx; 3195 b->next = NULL; 3196 /* Append to linked list so that we can preserve the order */ 3197 if (!head) ds->boundary = b; 3198 while (head) { 3199 if (!head->next) { 3200 head->next = b; 3201 head = b; 3202 } 3203 head = head->next; 3204 ++n; 3205 } 3206 if (bd) {PetscValidIntPointer(bd, 12); *bd = n;} 3207 PetscFunctionReturn(0); 3208 } 3209 3210 /*@C 3211 PetscDSAddBoundaryByName - Add a boundary condition to the model. The pointwise functions are used to provide boundary values for essential boundary conditions. In FEM, they are acting upon by dual basis functionals to generate FEM coefficients which are fixed. Natural boundary conditions signal to PETSc that boundary integrals should be performaed, using the kernels from PetscDSSetBdResidual(). 3212 3213 Collective on ds 3214 3215 Input Parameters: 3216 + ds - The PetscDS object 3217 . type - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann) 3218 . name - The BC name 3219 . lname - The naem of the label defining constrained points 3220 . Nv - The number of DMLabel values for constrained points 3221 . values - An array of label values for constrained points 3222 . field - The field to constrain 3223 . Nc - The number of constrained field components (0 will constrain all fields) 3224 . comps - An array of constrained component numbers 3225 . bcFunc - A pointwise function giving boundary values 3226 . bcFunc_t - A pointwise function giving the time derviative of the boundary values, or NULL 3227 - ctx - An optional user context for bcFunc 3228 3229 Output Parameters: 3230 - bd - The boundary number 3231 3232 Options Database Keys: 3233 + -bc_<boundary name> <num> - Overrides the boundary ids 3234 - -bc_<boundary name>_comp <num> - Overrides the boundary components 3235 3236 Note: 3237 This function should only be used with DMForest currently, since labels cannot be defined before the underlygin Plex is built. 3238 3239 Both bcFunc abd bcFunc_t will depend on the boundary condition type. If the type if DM_BC_ESSENTIAL, Then the calling sequence is: 3240 3241 $ bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[]) 3242 3243 If the type is DM_BC_ESSENTIAL_FIELD or other _FIELD value, then the calling sequence is: 3244 3245 $ bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, 3246 $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 3247 $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 3248 $ PetscReal time, const PetscReal x[], PetscScalar bcval[]) 3249 3250 + dim - the spatial dimension 3251 . Nf - the number of fields 3252 . uOff - the offset into u[] and u_t[] for each field 3253 . uOff_x - the offset into u_x[] for each field 3254 . u - each field evaluated at the current point 3255 . u_t - the time derivative of each field evaluated at the current point 3256 . u_x - the gradient of each field evaluated at the current point 3257 . aOff - the offset into a[] and a_t[] for each auxiliary field 3258 . aOff_x - the offset into a_x[] for each auxiliary field 3259 . a - each auxiliary field evaluated at the current point 3260 . a_t - the time derivative of each auxiliary field evaluated at the current point 3261 . a_x - the gradient of auxiliary each field evaluated at the current point 3262 . t - current time 3263 . x - coordinates of the current point 3264 . numConstants - number of constant parameters 3265 . constants - constant parameters 3266 - bcval - output values at the current point 3267 3268 Level: developer 3269 3270 .seealso: PetscDSAddBoundary(), PetscDSGetBoundary(), PetscDSSetResidual(), PetscDSSetBdResidual() 3271 @*/ 3272 PetscErrorCode PetscDSAddBoundaryByName(PetscDS ds, DMBoundaryConditionType type, const char name[], const char lname[], PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx, PetscInt *bd) 3273 { 3274 DSBoundary head = ds->boundary, b; 3275 PetscInt n = 0; 3276 PetscErrorCode ierr; 3277 3278 PetscFunctionBegin; 3279 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3280 PetscValidLogicalCollectiveEnum(ds, type, 2); 3281 PetscValidCharPointer(name, 3); 3282 PetscValidCharPointer(lname, 4); 3283 PetscValidLogicalCollectiveInt(ds, Nv, 5); 3284 PetscValidLogicalCollectiveInt(ds, field, 7); 3285 PetscValidLogicalCollectiveInt(ds, Nc, 8); 3286 ierr = PetscNew(&b);CHKERRQ(ierr); 3287 ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr); 3288 ierr = PetscWeakFormCreate(PETSC_COMM_SELF, &b->wf);CHKERRQ(ierr); 3289 ierr = PetscWeakFormSetNumFields(b->wf, ds->Nf);CHKERRQ(ierr); 3290 ierr = PetscMalloc1(Nv, &b->values);CHKERRQ(ierr); 3291 if (Nv) {ierr = PetscArraycpy(b->values, values, Nv);CHKERRQ(ierr);} 3292 ierr = PetscMalloc1(Nc, &b->comps);CHKERRQ(ierr); 3293 if (Nc) {ierr = PetscArraycpy(b->comps, comps, Nc);CHKERRQ(ierr);} 3294 ierr = PetscStrallocpy(lname, (char **) &b->lname);CHKERRQ(ierr); 3295 b->type = type; 3296 b->label = NULL; 3297 b->Nv = Nv; 3298 b->field = field; 3299 b->Nc = Nc; 3300 b->func = bcFunc; 3301 b->func_t = bcFunc_t; 3302 b->ctx = ctx; 3303 b->next = NULL; 3304 /* Append to linked list so that we can preserve the order */ 3305 if (!head) ds->boundary = b; 3306 while (head) { 3307 if (!head->next) { 3308 head->next = b; 3309 head = b; 3310 } 3311 head = head->next; 3312 ++n; 3313 } 3314 if (bd) {PetscValidIntPointer(bd, 12); *bd = n;} 3315 PetscFunctionReturn(0); 3316 } 3317 3318 /*@C 3319 PetscDSUpdateBoundary - Change a boundary condition for the model. The pointwise functions are used to provide boundary values for essential boundary conditions. In FEM, they are acting upon by dual basis functionals to generate FEM coefficients which are fixed. Natural boundary conditions signal to PETSc that boundary integrals should be performaed, using the kernels from PetscDSSetBdResidual(). 3320 3321 Input Parameters: 3322 + ds - The PetscDS object 3323 . bd - The boundary condition number 3324 . type - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann) 3325 . name - The BC name 3326 . label - The label defining constrained points 3327 . Nv - The number of DMLabel ids for constrained points 3328 . values - An array of ids for constrained points 3329 . field - The field to constrain 3330 . Nc - The number of constrained field components 3331 . comps - An array of constrained component numbers 3332 . bcFunc - A pointwise function giving boundary values 3333 . bcFunc_t - A pointwise function giving the time derviative of the boundary values, or NULL 3334 - ctx - An optional user context for bcFunc 3335 3336 Note: 3337 The boundary condition number is the order in which it was registered. The user can get the number of boundary conditions from PetscDSGetNumBoundary(). See PetscDSAddBoundary() for a description of the calling sequences for the callbacks. 3338 3339 Level: developer 3340 3341 .seealso: PetscDSAddBoundary(), PetscDSGetBoundary(), PetscDSGetNumBoundary() 3342 @*/ 3343 PetscErrorCode PetscDSUpdateBoundary(PetscDS ds, PetscInt bd, DMBoundaryConditionType type, const char name[], DMLabel label, PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx) 3344 { 3345 DSBoundary b = ds->boundary; 3346 PetscInt n = 0; 3347 PetscErrorCode ierr; 3348 3349 PetscFunctionBegin; 3350 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3351 while (b) { 3352 if (n == bd) break; 3353 b = b->next; 3354 ++n; 3355 } 3356 if (!b) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n); 3357 if (name) { 3358 ierr = PetscFree(b->name);CHKERRQ(ierr); 3359 ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr); 3360 } 3361 b->type = type; 3362 if (label) { 3363 const char *name; 3364 3365 b->label = label; 3366 ierr = PetscFree(b->lname);CHKERRQ(ierr); 3367 ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr); 3368 ierr = PetscStrallocpy(name, (char **) &b->lname);CHKERRQ(ierr); 3369 } 3370 if (Nv >= 0) { 3371 b->Nv = Nv; 3372 ierr = PetscFree(b->values);CHKERRQ(ierr); 3373 ierr = PetscMalloc1(Nv, &b->values);CHKERRQ(ierr); 3374 if (Nv) {ierr = PetscArraycpy(b->values, values, Nv);CHKERRQ(ierr);} 3375 } 3376 if (field >= 0) b->field = field; 3377 if (Nc >= 0) { 3378 b->Nc = Nc; 3379 ierr = PetscFree(b->comps);CHKERRQ(ierr); 3380 ierr = PetscMalloc1(Nc, &b->comps);CHKERRQ(ierr); 3381 if (Nc) {ierr = PetscArraycpy(b->comps, comps, Nc);CHKERRQ(ierr);} 3382 } 3383 if (bcFunc) b->func = bcFunc; 3384 if (bcFunc_t) b->func_t = bcFunc_t; 3385 if (ctx) b->ctx = ctx; 3386 PetscFunctionReturn(0); 3387 } 3388 3389 /*@ 3390 PetscDSGetNumBoundary - Get the number of registered BC 3391 3392 Input Parameters: 3393 . ds - The PetscDS object 3394 3395 Output Parameters: 3396 . numBd - The number of BC 3397 3398 Level: intermediate 3399 3400 .seealso: PetscDSAddBoundary(), PetscDSGetBoundary() 3401 @*/ 3402 PetscErrorCode PetscDSGetNumBoundary(PetscDS ds, PetscInt *numBd) 3403 { 3404 DSBoundary b = ds->boundary; 3405 3406 PetscFunctionBegin; 3407 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3408 PetscValidPointer(numBd, 2); 3409 *numBd = 0; 3410 while (b) {++(*numBd); b = b->next;} 3411 PetscFunctionReturn(0); 3412 } 3413 3414 /*@C 3415 PetscDSGetBoundary - Gets a boundary condition to the model 3416 3417 Input Parameters: 3418 + ds - The PetscDS object 3419 - bd - The BC number 3420 3421 Output Parameters: 3422 + wf - The PetscWeakForm holding the pointwise functions 3423 . type - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann) 3424 . name - The BC name 3425 . label - The label defining constrained points 3426 . Nv - The number of DMLabel ids for constrained points 3427 . values - An array of ids for constrained points 3428 . field - The field to constrain 3429 . Nc - The number of constrained field components 3430 . comps - An array of constrained component numbers 3431 . bcFunc - A pointwise function giving boundary values 3432 . bcFunc_t - A pointwise function giving the time derviative of the boundary values 3433 - ctx - An optional user context for bcFunc 3434 3435 Options Database Keys: 3436 + -bc_<boundary name> <num> - Overrides the boundary ids 3437 - -bc_<boundary name>_comp <num> - Overrides the boundary components 3438 3439 Level: developer 3440 3441 .seealso: PetscDSAddBoundary() 3442 @*/ 3443 PetscErrorCode PetscDSGetBoundary(PetscDS ds, PetscInt bd, PetscWeakForm *wf, DMBoundaryConditionType *type, const char *name[], DMLabel *label, PetscInt *Nv, const PetscInt *values[], PetscInt *field, PetscInt *Nc, const PetscInt *comps[], void (**func)(void), void (**func_t)(void), void **ctx) 3444 { 3445 DSBoundary b = ds->boundary; 3446 PetscInt n = 0; 3447 3448 PetscFunctionBegin; 3449 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3450 while (b) { 3451 if (n == bd) break; 3452 b = b->next; 3453 ++n; 3454 } 3455 if (!b) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n); 3456 if (wf) { 3457 PetscValidPointer(wf, 3); 3458 *wf = b->wf; 3459 } 3460 if (type) { 3461 PetscValidPointer(type, 4); 3462 *type = b->type; 3463 } 3464 if (name) { 3465 PetscValidPointer(name, 5); 3466 *name = b->name; 3467 } 3468 if (label) { 3469 PetscValidPointer(label, 6); 3470 *label = b->label; 3471 } 3472 if (Nv) { 3473 PetscValidIntPointer(Nv, 7); 3474 *Nv = b->Nv; 3475 } 3476 if (values) { 3477 PetscValidPointer(values, 8); 3478 *values = b->values; 3479 } 3480 if (field) { 3481 PetscValidIntPointer(field, 9); 3482 *field = b->field; 3483 } 3484 if (Nc) { 3485 PetscValidIntPointer(Nc, 10); 3486 *Nc = b->Nc; 3487 } 3488 if (comps) { 3489 PetscValidPointer(comps, 11); 3490 *comps = b->comps; 3491 } 3492 if (func) { 3493 PetscValidPointer(func, 12); 3494 *func = b->func; 3495 } 3496 if (func_t) { 3497 PetscValidPointer(func_t, 13); 3498 *func_t = b->func_t; 3499 } 3500 if (ctx) { 3501 PetscValidPointer(ctx, 14); 3502 *ctx = b->ctx; 3503 } 3504 PetscFunctionReturn(0); 3505 } 3506 3507 static PetscErrorCode DSBoundaryDuplicate_Internal(DSBoundary b, DSBoundary *bNew) 3508 { 3509 PetscErrorCode ierr; 3510 3511 PetscFunctionBegin; 3512 ierr = PetscNew(bNew);CHKERRQ(ierr); 3513 ierr = PetscWeakFormCreate(PETSC_COMM_SELF, &(*bNew)->wf);CHKERRQ(ierr); 3514 ierr = PetscWeakFormCopy(b->wf, (*bNew)->wf);CHKERRQ(ierr); 3515 ierr = PetscStrallocpy(b->name,(char **) &((*bNew)->name));CHKERRQ(ierr); 3516 ierr = PetscStrallocpy(b->lname,(char **) &((*bNew)->lname));CHKERRQ(ierr); 3517 (*bNew)->type = b->type; 3518 (*bNew)->label = b->label; 3519 (*bNew)->Nv = b->Nv; 3520 ierr = PetscMalloc1(b->Nv, &(*bNew)->values);CHKERRQ(ierr); 3521 ierr = PetscArraycpy((*bNew)->values, b->values, b->Nv);CHKERRQ(ierr); 3522 (*bNew)->field = b->field; 3523 (*bNew)->Nc = b->Nc; 3524 ierr = PetscMalloc1(b->Nc, &(*bNew)->comps);CHKERRQ(ierr); 3525 ierr = PetscArraycpy((*bNew)->comps, b->comps, b->Nc);CHKERRQ(ierr); 3526 (*bNew)->func = b->func; 3527 (*bNew)->func_t = b->func_t; 3528 (*bNew)->ctx = b->ctx; 3529 PetscFunctionReturn(0); 3530 } 3531 3532 /*@ 3533 PetscDSCopyBoundary - Copy all boundary condition objects to the new problem 3534 3535 Not collective 3536 3537 Input Parameters: 3538 + ds - The source PetscDS object 3539 . numFields - The number of selected fields, or PETSC_DEFAULT for all fields 3540 - fields - The selected fields, or NULL for all fields 3541 3542 Output Parameter: 3543 . newds - The target PetscDS, now with a copy of the boundary conditions 3544 3545 Level: intermediate 3546 3547 .seealso: PetscDSCopyEquations(), PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate() 3548 @*/ 3549 PetscErrorCode PetscDSCopyBoundary(PetscDS ds, PetscInt numFields, const PetscInt fields[], PetscDS newds) 3550 { 3551 DSBoundary b, *lastnext; 3552 PetscErrorCode ierr; 3553 3554 PetscFunctionBegin; 3555 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3556 PetscValidHeaderSpecific(newds, PETSCDS_CLASSID, 4); 3557 if (ds == newds) PetscFunctionReturn(0); 3558 ierr = PetscDSDestroyBoundary(newds);CHKERRQ(ierr); 3559 lastnext = &(newds->boundary); 3560 for (b = ds->boundary; b; b = b->next) { 3561 DSBoundary bNew; 3562 PetscInt fieldNew = -1; 3563 3564 if (numFields > 0 && fields) { 3565 PetscInt f; 3566 3567 for (f = 0; f < numFields; ++f) if (b->field == fields[f]) break; 3568 if (f == numFields) continue; 3569 fieldNew = f; 3570 } 3571 ierr = DSBoundaryDuplicate_Internal(b, &bNew);CHKERRQ(ierr); 3572 bNew->field = fieldNew < 0 ? b->field : fieldNew; 3573 *lastnext = bNew; 3574 lastnext = &(bNew->next); 3575 } 3576 PetscFunctionReturn(0); 3577 } 3578 3579 /*@ 3580 PetscDSDestroyBoundary - Remove all DMBoundary objects from the PetscDS 3581 3582 Not collective 3583 3584 Input Parameter: 3585 . ds - The PetscDS object 3586 3587 Level: intermediate 3588 3589 .seealso: PetscDSCopyBoundary(), PetscDSCopyEquations() 3590 @*/ 3591 PetscErrorCode PetscDSDestroyBoundary(PetscDS ds) 3592 { 3593 DSBoundary next = ds->boundary; 3594 PetscErrorCode ierr; 3595 3596 PetscFunctionBegin; 3597 while (next) { 3598 DSBoundary b = next; 3599 3600 next = b->next; 3601 ierr = PetscWeakFormDestroy(&b->wf);CHKERRQ(ierr); 3602 ierr = PetscFree(b->name);CHKERRQ(ierr); 3603 ierr = PetscFree(b->lname);CHKERRQ(ierr); 3604 ierr = PetscFree(b->values);CHKERRQ(ierr); 3605 ierr = PetscFree(b->comps);CHKERRQ(ierr); 3606 ierr = PetscFree(b);CHKERRQ(ierr); 3607 } 3608 PetscFunctionReturn(0); 3609 } 3610 3611 /*@ 3612 PetscDSSelectDiscretizations - Copy discretizations to the new problem with different field layout 3613 3614 Not collective 3615 3616 Input Parameter: 3617 + prob - The PetscDS object 3618 . numFields - Number of new fields 3619 - fields - Old field number for each new field 3620 3621 Output Parameter: 3622 . newprob - The PetscDS copy 3623 3624 Level: intermediate 3625 3626 .seealso: PetscDSSelectEquations(), PetscDSCopyBoundary(), PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate() 3627 @*/ 3628 PetscErrorCode PetscDSSelectDiscretizations(PetscDS prob, PetscInt numFields, const PetscInt fields[], PetscDS newprob) 3629 { 3630 PetscInt Nf, Nfn, fn; 3631 PetscErrorCode ierr; 3632 3633 PetscFunctionBegin; 3634 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3635 if (fields) PetscValidPointer(fields, 3); 3636 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 4); 3637 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 3638 ierr = PetscDSGetNumFields(newprob, &Nfn);CHKERRQ(ierr); 3639 numFields = numFields < 0 ? Nf : numFields; 3640 for (fn = 0; fn < numFields; ++fn) { 3641 const PetscInt f = fields ? fields[fn] : fn; 3642 PetscObject disc; 3643 3644 if (f >= Nf) continue; 3645 ierr = PetscDSGetDiscretization(prob, f, &disc);CHKERRQ(ierr); 3646 ierr = PetscDSSetDiscretization(newprob, fn, disc);CHKERRQ(ierr); 3647 } 3648 PetscFunctionReturn(0); 3649 } 3650 3651 /*@ 3652 PetscDSSelectEquations - Copy pointwise function pointers to the new problem with different field layout 3653 3654 Not collective 3655 3656 Input Parameter: 3657 + prob - The PetscDS object 3658 . numFields - Number of new fields 3659 - fields - Old field number for each new field 3660 3661 Output Parameter: 3662 . newprob - The PetscDS copy 3663 3664 Level: intermediate 3665 3666 .seealso: PetscDSSelectDiscretizations(), PetscDSCopyBoundary(), PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate() 3667 @*/ 3668 PetscErrorCode PetscDSSelectEquations(PetscDS prob, PetscInt numFields, const PetscInt fields[], PetscDS newprob) 3669 { 3670 PetscInt Nf, Nfn, fn, gn; 3671 PetscErrorCode ierr; 3672 3673 PetscFunctionBegin; 3674 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3675 if (fields) PetscValidPointer(fields, 3); 3676 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 4); 3677 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 3678 ierr = PetscDSGetNumFields(newprob, &Nfn);CHKERRQ(ierr); 3679 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); 3680 for (fn = 0; fn < numFields; ++fn) { 3681 const PetscInt f = fields ? fields[fn] : fn; 3682 PetscPointFunc obj; 3683 PetscPointFunc f0, f1; 3684 PetscBdPointFunc f0Bd, f1Bd; 3685 PetscRiemannFunc r; 3686 3687 if (f >= Nf) continue; 3688 ierr = PetscDSGetObjective(prob, f, &obj);CHKERRQ(ierr); 3689 ierr = PetscDSGetResidual(prob, f, &f0, &f1);CHKERRQ(ierr); 3690 ierr = PetscDSGetBdResidual(prob, f, &f0Bd, &f1Bd);CHKERRQ(ierr); 3691 ierr = PetscDSGetRiemannSolver(prob, f, &r);CHKERRQ(ierr); 3692 ierr = PetscDSSetObjective(newprob, fn, obj);CHKERRQ(ierr); 3693 ierr = PetscDSSetResidual(newprob, fn, f0, f1);CHKERRQ(ierr); 3694 ierr = PetscDSSetBdResidual(newprob, fn, f0Bd, f1Bd);CHKERRQ(ierr); 3695 ierr = PetscDSSetRiemannSolver(newprob, fn, r);CHKERRQ(ierr); 3696 for (gn = 0; gn < numFields; ++gn) { 3697 const PetscInt g = fields ? fields[gn] : gn; 3698 PetscPointJac g0, g1, g2, g3; 3699 PetscPointJac g0p, g1p, g2p, g3p; 3700 PetscBdPointJac g0Bd, g1Bd, g2Bd, g3Bd; 3701 3702 if (g >= Nf) continue; 3703 ierr = PetscDSGetJacobian(prob, f, g, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 3704 ierr = PetscDSGetJacobianPreconditioner(prob, f, g, &g0p, &g1p, &g2p, &g3p);CHKERRQ(ierr); 3705 ierr = PetscDSGetBdJacobian(prob, f, g, &g0Bd, &g1Bd, &g2Bd, &g3Bd);CHKERRQ(ierr); 3706 ierr = PetscDSSetJacobian(newprob, fn, gn, g0, g1, g2, g3);CHKERRQ(ierr); 3707 ierr = PetscDSSetJacobianPreconditioner(newprob, fn, gn, g0p, g1p, g2p, g3p);CHKERRQ(ierr); 3708 ierr = PetscDSSetBdJacobian(newprob, fn, gn, g0Bd, g1Bd, g2Bd, g3Bd);CHKERRQ(ierr); 3709 } 3710 } 3711 PetscFunctionReturn(0); 3712 } 3713 3714 /*@ 3715 PetscDSCopyEquations - Copy all pointwise function pointers to the new problem 3716 3717 Not collective 3718 3719 Input Parameter: 3720 . prob - The PetscDS object 3721 3722 Output Parameter: 3723 . newprob - The PetscDS copy 3724 3725 Level: intermediate 3726 3727 .seealso: PetscDSCopyBoundary(), PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate() 3728 @*/ 3729 PetscErrorCode PetscDSCopyEquations(PetscDS prob, PetscDS newprob) 3730 { 3731 PetscInt Nf, Ng; 3732 PetscErrorCode ierr; 3733 3734 PetscFunctionBegin; 3735 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3736 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2); 3737 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 3738 ierr = PetscDSGetNumFields(newprob, &Ng);CHKERRQ(ierr); 3739 if (Nf != Ng) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_SIZ, "Number of fields must match %D != %D", Nf, Ng); 3740 ierr = PetscDSSelectEquations(prob, Nf, NULL, newprob);CHKERRQ(ierr); 3741 PetscFunctionReturn(0); 3742 } 3743 3744 /*@ 3745 PetscDSCopyConstants - Copy all constants to the new problem 3746 3747 Not collective 3748 3749 Input Parameter: 3750 . prob - The PetscDS object 3751 3752 Output Parameter: 3753 . newprob - The PetscDS copy 3754 3755 Level: intermediate 3756 3757 .seealso: PetscDSCopyBoundary(), PetscDSCopyEquations(), PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate() 3758 @*/ 3759 PetscErrorCode PetscDSCopyConstants(PetscDS prob, PetscDS newprob) 3760 { 3761 PetscInt Nc; 3762 const PetscScalar *constants; 3763 PetscErrorCode ierr; 3764 3765 PetscFunctionBegin; 3766 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3767 PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2); 3768 ierr = PetscDSGetConstants(prob, &Nc, &constants);CHKERRQ(ierr); 3769 ierr = PetscDSSetConstants(newprob, Nc, (PetscScalar *) constants);CHKERRQ(ierr); 3770 PetscFunctionReturn(0); 3771 } 3772 3773 /*@ 3774 PetscDSCopyExactSolutions - Copy all exact solutions to the new problem 3775 3776 Not collective 3777 3778 Input Parameter: 3779 . ds - The PetscDS object 3780 3781 Output Parameter: 3782 . newds - The PetscDS copy 3783 3784 Level: intermediate 3785 3786 .seealso: PetscDSCopyBoundary(), PetscDSCopyEquations(), PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate() 3787 @*/ 3788 PetscErrorCode PetscDSCopyExactSolutions(PetscDS ds, PetscDS newds) 3789 { 3790 PetscSimplePointFunc sol; 3791 void *ctx; 3792 PetscInt Nf, f; 3793 PetscErrorCode ierr; 3794 3795 PetscFunctionBegin; 3796 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3797 PetscValidHeaderSpecific(newds, PETSCDS_CLASSID, 2); 3798 ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 3799 for (f = 0; f < Nf; ++f) { 3800 ierr = PetscDSGetExactSolution(ds, f, &sol, &ctx);CHKERRQ(ierr); 3801 ierr = PetscDSSetExactSolution(newds, f, sol, ctx);CHKERRQ(ierr); 3802 ierr = PetscDSGetExactSolutionTimeDerivative(ds, f, &sol, &ctx);CHKERRQ(ierr); 3803 ierr = PetscDSSetExactSolutionTimeDerivative(newds, f, sol, ctx);CHKERRQ(ierr); 3804 } 3805 PetscFunctionReturn(0); 3806 } 3807 3808 PetscErrorCode PetscDSGetHeightSubspace(PetscDS prob, PetscInt height, PetscDS *subprob) 3809 { 3810 PetscInt dim, Nf, f; 3811 PetscErrorCode ierr; 3812 3813 PetscFunctionBegin; 3814 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 3815 PetscValidPointer(subprob, 3); 3816 if (height == 0) {*subprob = prob; PetscFunctionReturn(0);} 3817 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 3818 ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr); 3819 if (height > dim) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_OUTOFRANGE, "DS can only handle height in [0, %D], not %D", dim, height); 3820 if (!prob->subprobs) {ierr = PetscCalloc1(dim, &prob->subprobs);CHKERRQ(ierr);} 3821 if (!prob->subprobs[height-1]) { 3822 PetscInt cdim; 3823 3824 ierr = PetscDSCreate(PetscObjectComm((PetscObject) prob), &prob->subprobs[height-1]);CHKERRQ(ierr); 3825 ierr = PetscDSGetCoordinateDimension(prob, &cdim);CHKERRQ(ierr); 3826 ierr = PetscDSSetCoordinateDimension(prob->subprobs[height-1], cdim);CHKERRQ(ierr); 3827 for (f = 0; f < Nf; ++f) { 3828 PetscFE subfe; 3829 PetscObject obj; 3830 PetscClassId id; 3831 3832 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 3833 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 3834 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetHeightSubspace((PetscFE) obj, height, &subfe);CHKERRQ(ierr);} 3835 else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported discretization type for field %d", f); 3836 ierr = PetscDSSetDiscretization(prob->subprobs[height-1], f, (PetscObject) subfe);CHKERRQ(ierr); 3837 } 3838 } 3839 *subprob = prob->subprobs[height-1]; 3840 PetscFunctionReturn(0); 3841 } 3842 3843 PetscErrorCode PetscDSGetDiscType_Internal(PetscDS ds, PetscInt f, PetscDiscType *disctype) 3844 { 3845 PetscObject obj; 3846 PetscClassId id; 3847 PetscInt Nf; 3848 PetscErrorCode ierr; 3849 3850 PetscFunctionBegin; 3851 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3852 PetscValidPointer(disctype, 3); 3853 *disctype = PETSC_DISC_NONE; 3854 ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 3855 if (f >= Nf) SETERRQ2(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_SIZ, "Field %D must be in [0, %D)", f, Nf); 3856 ierr = PetscDSGetDiscretization(ds, f, &obj);CHKERRQ(ierr); 3857 if (obj) { 3858 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 3859 if (id == PETSCFE_CLASSID) *disctype = PETSC_DISC_FE; 3860 else *disctype = PETSC_DISC_FV; 3861 } 3862 PetscFunctionReturn(0); 3863 } 3864 3865 static PetscErrorCode PetscDSDestroy_Basic(PetscDS ds) 3866 { 3867 PetscErrorCode ierr; 3868 3869 PetscFunctionBegin; 3870 ierr = PetscFree(ds->data);CHKERRQ(ierr); 3871 PetscFunctionReturn(0); 3872 } 3873 3874 static PetscErrorCode PetscDSInitialize_Basic(PetscDS ds) 3875 { 3876 PetscFunctionBegin; 3877 ds->ops->setfromoptions = NULL; 3878 ds->ops->setup = NULL; 3879 ds->ops->view = NULL; 3880 ds->ops->destroy = PetscDSDestroy_Basic; 3881 PetscFunctionReturn(0); 3882 } 3883 3884 /*MC 3885 PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions 3886 3887 Level: intermediate 3888 3889 .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType() 3890 M*/ 3891 3892 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS ds) 3893 { 3894 PetscDS_Basic *b; 3895 PetscErrorCode ierr; 3896 3897 PetscFunctionBegin; 3898 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 3899 ierr = PetscNewLog(ds, &b);CHKERRQ(ierr); 3900 ds->data = b; 3901 3902 ierr = PetscDSInitialize_Basic(ds);CHKERRQ(ierr); 3903 PetscFunctionReturn(0); 3904 } 3905