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