1 #include <petsc/private/petscdsimpl.h> /*I "petscds.h" I*/ 2 3 PetscClassId PETSCWEAKFORM_CLASSID = 0; 4 5 const char *const PetscWeakFormKinds[] = {"objective", "residual_f0", "residual_f1", "jacobian_g0", "jacobian_g1", "jacobian_g2", "jacobian_g3", "jacobian_preconditioner_g0", "jacobian_preconditioner_g1", "jacobian_preconditioner_g2", "jacobian_preconditioner_g3", "dynamic_jacobian_g0", "dynamic_jacobian_g1", "dynamic_jacobian_g2", "dynamic_jacobian_g3", "boundary_residual_f0", "boundary_residual_f1", "boundary_jacobian_g0", "boundary_jacobian_g1", "boundary_jacobian_g2", "boundary_jacobian_g3", "boundary_jacobian_preconditioner_g0", "boundary_jacobian_preconditioner_g1", "boundary_jacobian_preconditioner_g2", "boundary_jacobian_preconditioner_g3", "riemann_solver", "PetscWeakFormKind", "PETSC_WF_", NULL}; 6 7 static PetscErrorCode PetscChunkBufferCreate(size_t unitbytes, size_t expected, PetscChunkBuffer **buffer) 8 { 9 PetscFunctionBegin; 10 PetscCall(PetscNew(buffer)); 11 PetscCall(PetscCalloc1(expected * unitbytes, &(*buffer)->array)); 12 (*buffer)->size = expected; 13 (*buffer)->unitbytes = unitbytes; 14 (*buffer)->alloc = expected * unitbytes; 15 PetscFunctionReturn(0); 16 } 17 18 static PetscErrorCode PetscChunkBufferDuplicate(PetscChunkBuffer *buffer, PetscChunkBuffer **bufferNew) 19 { 20 PetscFunctionBegin; 21 PetscCall(PetscNew(bufferNew)); 22 PetscCall(PetscCalloc1(buffer->size * buffer->unitbytes, &(*bufferNew)->array)); 23 PetscCall(PetscMemcpy((*bufferNew)->array, buffer->array, buffer->size * buffer->unitbytes)); 24 (*bufferNew)->size = buffer->size; 25 (*bufferNew)->unitbytes = buffer->unitbytes; 26 (*bufferNew)->alloc = buffer->size * buffer->unitbytes; 27 PetscFunctionReturn(0); 28 } 29 30 static PetscErrorCode PetscChunkBufferDestroy(PetscChunkBuffer **buffer) 31 { 32 PetscFunctionBegin; 33 PetscCall(PetscFree((*buffer)->array)); 34 PetscCall(PetscFree(*buffer)); 35 PetscFunctionReturn(0); 36 } 37 38 static PetscErrorCode PetscChunkBufferCreateChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk) 39 { 40 PetscFunctionBegin; 41 if ((buffer->size + size) * buffer->unitbytes > buffer->alloc) { 42 char *tmp; 43 44 if (!buffer->alloc) buffer->alloc = (buffer->size + size) * buffer->unitbytes; 45 while ((buffer->size + size) * buffer->unitbytes > buffer->alloc) buffer->alloc *= 2; 46 PetscCall(PetscMalloc(buffer->alloc, &tmp)); 47 PetscCall(PetscMemcpy(tmp, buffer->array, buffer->size * buffer->unitbytes)); 48 PetscCall(PetscFree(buffer->array)); 49 buffer->array = tmp; 50 } 51 chunk->start = buffer->size * buffer->unitbytes; 52 chunk->size = size; 53 chunk->reserved = size; 54 buffer->size += size; 55 PetscFunctionReturn(0); 56 } 57 58 static PetscErrorCode PetscChunkBufferEnlargeChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk) 59 { 60 size_t siz = size; 61 62 PetscFunctionBegin; 63 if (chunk->size + size > chunk->reserved) { 64 PetscChunk newchunk; 65 PetscInt reserved = chunk->size; 66 67 /* TODO Here if we had a chunk list, we could update them all to reclaim unused space */ 68 while (reserved < chunk->size + size) reserved *= 2; 69 PetscCall(PetscChunkBufferCreateChunk(buffer, (size_t)reserved, &newchunk)); 70 newchunk.size = chunk->size + size; 71 PetscCall(PetscMemcpy(&buffer->array[newchunk.start], &buffer->array[chunk->start], chunk->size * buffer->unitbytes)); 72 *chunk = newchunk; 73 } else { 74 chunk->size += siz; 75 } 76 PetscFunctionReturn(0); 77 } 78 79 /*@C 80 PetscFormKeySort - Sorts an array of `PetscFormKey` in place in increasing order. 81 82 Not Collective 83 84 Input Parameters: 85 + n - number of values 86 - X - array of PetscFormKey 87 88 Level: intermediate 89 90 .seealso: `PetscFormKey`, `PetscIntSortSemiOrdered()`, `PetscSortInt()` 91 @*/ 92 PetscErrorCode PetscFormKeySort(PetscInt n, PetscFormKey arr[]) 93 { 94 PetscFunctionBegin; 95 if (n <= 1) PetscFunctionReturn(0); 96 PetscValidPointer(arr, 2); 97 PetscCall(PetscTimSort(n, arr, sizeof(PetscFormKey), Compare_PetscFormKey_Private, NULL)); 98 PetscFunctionReturn(0); 99 } 100 101 PetscErrorCode PetscWeakFormGetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt *n, void (***func)(void)) 102 { 103 PetscFormKey key; 104 PetscChunk chunk; 105 106 PetscFunctionBegin; 107 key.label = label; 108 key.value = value; 109 key.field = f; 110 key.part = part; 111 PetscCall(PetscHMapFormGet(ht, key, &chunk)); 112 if (chunk.size < 0) { 113 *n = 0; 114 *func = NULL; 115 } else { 116 *n = chunk.size; 117 *func = (void (**)(void)) & wf->funcs->array[chunk.start]; 118 } 119 PetscFunctionReturn(0); 120 } 121 122 /* A NULL argument for func causes this to clear the key */ 123 PetscErrorCode PetscWeakFormSetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt n, void (**func)(void)) 124 { 125 PetscFormKey key; 126 PetscChunk chunk; 127 PetscInt i; 128 129 PetscFunctionBegin; 130 key.label = label; 131 key.value = value; 132 key.field = f; 133 key.part = part; 134 if (!func) { 135 PetscCall(PetscHMapFormDel(ht, key)); 136 PetscFunctionReturn(0); 137 } else PetscCall(PetscHMapFormGet(ht, key, &chunk)); 138 if (chunk.size < 0) { 139 PetscCall(PetscChunkBufferCreateChunk(wf->funcs, n, &chunk)); 140 PetscCall(PetscHMapFormSet(ht, key, chunk)); 141 } else if (chunk.size <= n) { 142 PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, n - chunk.size, &chunk)); 143 PetscCall(PetscHMapFormSet(ht, key, chunk)); 144 } 145 for (i = 0; i < n; ++i) ((void (**)(void)) & wf->funcs->array[chunk.start])[i] = func[i]; 146 PetscFunctionReturn(0); 147 } 148 149 PetscErrorCode PetscWeakFormAddFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, void (*func)(void)) 150 { 151 PetscFormKey key; 152 PetscChunk chunk; 153 154 PetscFunctionBegin; 155 if (!func) PetscFunctionReturn(0); 156 key.label = label; 157 key.value = value; 158 key.field = f; 159 key.part = part; 160 PetscCall(PetscHMapFormGet(ht, key, &chunk)); 161 if (chunk.size < 0) { 162 PetscCall(PetscChunkBufferCreateChunk(wf->funcs, 1, &chunk)); 163 PetscCall(PetscHMapFormSet(ht, key, chunk)); 164 ((void (**)(void)) & wf->funcs->array[chunk.start])[0] = func; 165 } else { 166 PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, 1, &chunk)); 167 PetscCall(PetscHMapFormSet(ht, key, chunk)); 168 ((void (**)(void)) & wf->funcs->array[chunk.start])[chunk.size - 1] = func; 169 } 170 PetscFunctionReturn(0); 171 } 172 173 PetscErrorCode PetscWeakFormGetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, void (**func)(void)) 174 { 175 PetscFormKey key; 176 PetscChunk chunk; 177 178 PetscFunctionBegin; 179 key.label = label; 180 key.value = value; 181 key.field = f; 182 key.part = part; 183 PetscCall(PetscHMapFormGet(ht, key, &chunk)); 184 if (chunk.size < 0) { 185 *func = NULL; 186 } else { 187 PetscCheck(ind < chunk.size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", ind, chunk.size); 188 *func = ((void (**)(void)) & wf->funcs->array[chunk.start])[ind]; 189 } 190 PetscFunctionReturn(0); 191 } 192 193 /* Ignore a NULL func */ 194 PetscErrorCode PetscWeakFormSetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, void (*func)(void)) 195 { 196 PetscFormKey key; 197 PetscChunk chunk; 198 199 PetscFunctionBegin; 200 if (!func) PetscFunctionReturn(0); 201 key.label = label; 202 key.value = value; 203 key.field = f; 204 key.part = part; 205 PetscCall(PetscHMapFormGet(ht, key, &chunk)); 206 if (chunk.size < 0) { 207 PetscCall(PetscChunkBufferCreateChunk(wf->funcs, ind + 1, &chunk)); 208 PetscCall(PetscHMapFormSet(ht, key, chunk)); 209 } else if (chunk.size <= ind) { 210 PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, ind - chunk.size + 1, &chunk)); 211 PetscCall(PetscHMapFormSet(ht, key, chunk)); 212 } 213 ((void (**)(void)) & wf->funcs->array[chunk.start])[ind] = func; 214 PetscFunctionReturn(0); 215 } 216 217 PetscErrorCode PetscWeakFormClearIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind) 218 { 219 PetscFormKey key; 220 PetscChunk chunk; 221 222 PetscFunctionBegin; 223 key.label = label; 224 key.value = value; 225 key.field = f; 226 key.part = part; 227 PetscCall(PetscHMapFormGet(ht, key, &chunk)); 228 if (chunk.size < 0) { 229 PetscFunctionReturn(0); 230 } else if (!ind && chunk.size == 1) { 231 PetscCall(PetscHMapFormDel(ht, key)); 232 PetscFunctionReturn(0); 233 } else if (chunk.size <= ind) { 234 PetscFunctionReturn(0); 235 } 236 ((void (**)(void)) & wf->funcs->array[chunk.start])[ind] = NULL; 237 PetscFunctionReturn(0); 238 } 239 240 /*@ 241 PetscWeakFormCopy - Copy the pointwise functions to another `PetscWeakForm` 242 243 Not Collective 244 245 Input Parameter: 246 . wf - The original `PetscWeakForm` 247 248 Output Parameter: 249 . wfNew - The copy of the `PetscWeakForm` 250 251 Level: intermediate 252 253 .seealso: `PetscWeakForm`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 254 @*/ 255 PetscErrorCode PetscWeakFormCopy(PetscWeakForm wf, PetscWeakForm wfNew) 256 { 257 PetscInt f; 258 259 PetscFunctionBegin; 260 wfNew->Nf = wf->Nf; 261 PetscCall(PetscChunkBufferDestroy(&wfNew->funcs)); 262 PetscCall(PetscChunkBufferDuplicate(wf->funcs, &wfNew->funcs)); 263 for (f = 0; f < PETSC_NUM_WF; ++f) { 264 PetscCall(PetscHMapFormDestroy(&wfNew->form[f])); 265 PetscCall(PetscHMapFormDuplicate(wf->form[f], &wfNew->form[f])); 266 } 267 PetscFunctionReturn(0); 268 } 269 270 /*@ 271 PetscWeakFormClear - Clear all functions from the `PetscWeakForm` 272 273 Not Collective 274 275 Input Parameter: 276 . wf - The original `PetscWeakForm` 277 278 Level: intermediate 279 280 .seealso: `PetscWeakForm`, `PetscWeakFormCopy()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 281 @*/ 282 PetscErrorCode PetscWeakFormClear(PetscWeakForm wf) 283 { 284 PetscInt f; 285 286 PetscFunctionBegin; 287 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormClear(wf->form[f])); 288 PetscFunctionReturn(0); 289 } 290 291 static PetscErrorCode PetscWeakFormRewriteKeys_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label, PetscInt Nv, const PetscInt values[]) 292 { 293 PetscFormKey *keys; 294 PetscInt n, i, v, off = 0; 295 296 PetscFunctionBegin; 297 PetscCall(PetscHMapFormGetSize(hmap, &n)); 298 PetscCall(PetscMalloc1(n, &keys)); 299 PetscCall(PetscHMapFormGetKeys(hmap, &off, keys)); 300 for (i = 0; i < n; ++i) { 301 if (keys[i].label == label) { 302 PetscBool clear = PETSC_TRUE; 303 void (**funcs)(void); 304 PetscInt Nf; 305 306 PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 307 for (v = 0; v < Nv; ++v) { 308 PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, values[v], keys[i].field, keys[i].part, Nf, funcs)); 309 if (values[v] == keys[i].value) clear = PETSC_FALSE; 310 } 311 if (clear) PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL)); 312 } 313 } 314 PetscCall(PetscFree(keys)); 315 PetscFunctionReturn(0); 316 } 317 318 /*@C 319 PetscWeakFormRewriteKeys - Change any key on the given label to use the new set of label values 320 321 Not Collective 322 323 Input Parameters: 324 + wf - The original `PetscWeakForm` 325 . label - The label to change keys for 326 . Nv - The number of new label values 327 - values - The set of new values to relabel keys with 328 329 Level: intermediate 330 331 Note: 332 This is used internally when boundary label values are specified from the command line. 333 334 .seealso: `PetscWeakForm`, `DMLabel`, `PetscWeakFormReplaceLabel()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 335 @*/ 336 PetscErrorCode PetscWeakFormRewriteKeys(PetscWeakForm wf, DMLabel label, PetscInt Nv, const PetscInt values[]) 337 { 338 PetscInt f; 339 340 PetscFunctionBegin; 341 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormRewriteKeys_Internal(wf, wf->form[f], label, Nv, values)); 342 PetscFunctionReturn(0); 343 } 344 345 static PetscErrorCode PetscWeakFormReplaceLabel_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label) 346 { 347 PetscFormKey *keys; 348 PetscInt n, i, off = 0, maxFuncs = 0; 349 void (**tmpf)(void); 350 const char *name = NULL; 351 352 PetscFunctionBegin; 353 if (label) PetscCall(PetscObjectGetName((PetscObject)label, &name)); 354 PetscCall(PetscHMapFormGetSize(hmap, &n)); 355 PetscCall(PetscMalloc1(n, &keys)); 356 PetscCall(PetscHMapFormGetKeys(hmap, &off, keys)); 357 for (i = 0; i < n; ++i) { 358 PetscBool match = PETSC_FALSE; 359 const char *lname = NULL; 360 361 if (label == keys[i].label) continue; 362 if (keys[i].label) PetscCall(PetscObjectGetName((PetscObject)keys[i].label, &lname)); 363 PetscCall(PetscStrcmp(name, lname, &match)); 364 if ((!name && !lname) || match) { 365 void (**funcs)(void); 366 PetscInt Nf; 367 368 PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 369 maxFuncs = PetscMax(maxFuncs, Nf); 370 } 371 } 372 /* Need temp space because chunk buffer can be reallocated in SetFunction() call */ 373 PetscCall(PetscMalloc1(maxFuncs, &tmpf)); 374 for (i = 0; i < n; ++i) { 375 PetscBool match = PETSC_FALSE; 376 const char *lname = NULL; 377 378 if (label == keys[i].label) continue; 379 if (keys[i].label) PetscCall(PetscObjectGetName((PetscObject)keys[i].label, &lname)); 380 PetscCall(PetscStrcmp(name, lname, &match)); 381 if ((!name && !lname) || match) { 382 void (**funcs)(void); 383 PetscInt Nf, j; 384 385 PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 386 for (j = 0; j < Nf; ++j) tmpf[j] = funcs[j]; 387 PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, label, keys[i].value, keys[i].field, keys[i].part, Nf, tmpf)); 388 PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL)); 389 } 390 } 391 PetscCall(PetscFree(tmpf)); 392 PetscCall(PetscFree(keys)); 393 PetscFunctionReturn(0); 394 } 395 396 /*@C 397 PetscWeakFormReplaceLabel - Change any key on a label of the same name to use the new label 398 399 Not Collective 400 401 Input Parameters: 402 + wf - The original `PetscWeakForm` 403 - label - The label to change keys for 404 405 Level: intermediate 406 407 Note: 408 This is used internally when meshes are modified 409 410 .seealso: `PetscWeakForm`, `DMLabel`, `PetscWeakFormRewriteKeys()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 411 @*/ 412 PetscErrorCode PetscWeakFormReplaceLabel(PetscWeakForm wf, DMLabel label) 413 { 414 PetscInt f; 415 416 PetscFunctionBegin; 417 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormReplaceLabel_Internal(wf, wf->form[f], label)); 418 PetscFunctionReturn(0); 419 } 420 421 PetscErrorCode PetscWeakFormClearIndex(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscWeakFormKind kind, PetscInt ind) 422 { 423 PetscFunctionBegin; 424 PetscCall(PetscWeakFormClearIndexFunction_Private(wf, wf->form[kind], label, val, f, part, ind)); 425 PetscFunctionReturn(0); 426 } 427 428 PetscErrorCode PetscWeakFormGetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n, void (***obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 429 { 430 PetscFunctionBegin; 431 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (***)(void))obj)); 432 PetscFunctionReturn(0); 433 } 434 435 PetscErrorCode PetscWeakFormSetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n, void (**obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 436 { 437 PetscFunctionBegin; 438 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (**)(void))obj)); 439 PetscFunctionReturn(0); 440 } 441 442 PetscErrorCode PetscWeakFormAddObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, void (*obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 443 { 444 PetscFunctionBegin; 445 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, (void (*)(void))obj)); 446 PetscFunctionReturn(0); 447 } 448 449 PetscErrorCode PetscWeakFormGetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt ind, void (**obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 450 { 451 PetscFunctionBegin; 452 PetscCall(PetscWeakFormGetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (void (**)(void))obj)); 453 PetscFunctionReturn(0); 454 } 455 456 PetscErrorCode PetscWeakFormSetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt ind, void (*obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 457 { 458 PetscFunctionBegin; 459 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (void (*)(void))obj)); 460 PetscFunctionReturn(0); 461 } 462 463 PetscErrorCode PetscWeakFormGetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n0, void (***f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 464 { 465 PetscFunctionBegin; 466 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (***)(void))f0)); 467 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (***)(void))f1)); 468 PetscFunctionReturn(0); 469 } 470 471 PetscErrorCode PetscWeakFormAddResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 472 { 473 PetscFunctionBegin; 474 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, (void (*)(void))f0)); 475 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, (void (*)(void))f1)); 476 PetscFunctionReturn(0); 477 } 478 479 PetscErrorCode PetscWeakFormSetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n0, void (**f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 480 { 481 PetscFunctionBegin; 482 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (**)(void))f0)); 483 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (**)(void))f1)); 484 PetscFunctionReturn(0); 485 } 486 487 PetscErrorCode PetscWeakFormSetIndexResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt i0, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 488 { 489 PetscFunctionBegin; 490 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, i0, (void (*)(void))f0)); 491 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, i1, (void (*)(void))f1)); 492 PetscFunctionReturn(0); 493 } 494 495 PetscErrorCode PetscWeakFormGetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n0, void (***f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 496 { 497 PetscFunctionBegin; 498 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (***)(void))f0)); 499 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (***)(void))f1)); 500 PetscFunctionReturn(0); 501 } 502 503 PetscErrorCode PetscWeakFormAddBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 504 { 505 PetscFunctionBegin; 506 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, (void (*)(void))f0)); 507 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, (void (*)(void))f1)); 508 PetscFunctionReturn(0); 509 } 510 511 PetscErrorCode PetscWeakFormSetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n0, void (**f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 512 { 513 PetscFunctionBegin; 514 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (**)(void))f0)); 515 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (**)(void))f1)); 516 PetscFunctionReturn(0); 517 } 518 519 PetscErrorCode PetscWeakFormSetIndexBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt i0, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 520 { 521 PetscFunctionBegin; 522 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, i0, (void (*)(void))f0)); 523 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, i1, (void (*)(void))f1)); 524 PetscFunctionReturn(0); 525 } 526 527 PetscErrorCode PetscWeakFormHasJacobian(PetscWeakForm wf, PetscBool *hasJac) 528 { 529 PetscInt n0, n1, n2, n3; 530 531 PetscFunctionBegin; 532 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 533 PetscValidBoolPointer(hasJac, 2); 534 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G0], &n0)); 535 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G1], &n1)); 536 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G2], &n2)); 537 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G3], &n3)); 538 *hasJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 539 PetscFunctionReturn(0); 540 } 541 542 PetscErrorCode PetscWeakFormGetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 543 { 544 PetscInt find = f * wf->Nf + g; 545 546 PetscFunctionBegin; 547 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (***)(void))g0)); 548 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (***)(void))g1)); 549 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (***)(void))g2)); 550 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (***)(void))g3)); 551 PetscFunctionReturn(0); 552 } 553 554 PetscErrorCode PetscWeakFormAddJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 555 { 556 PetscInt find = f * wf->Nf + g; 557 558 PetscFunctionBegin; 559 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, (void (*)(void))g0)); 560 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, (void (*)(void))g1)); 561 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, (void (*)(void))g2)); 562 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, (void (*)(void))g3)); 563 PetscFunctionReturn(0); 564 } 565 566 PetscErrorCode PetscWeakFormSetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 567 { 568 PetscInt find = f * wf->Nf + g; 569 570 PetscFunctionBegin; 571 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (**)(void))g0)); 572 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (**)(void))g1)); 573 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (**)(void))g2)); 574 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (**)(void))g3)); 575 PetscFunctionReturn(0); 576 } 577 578 PetscErrorCode PetscWeakFormSetIndexJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 579 { 580 PetscInt find = f * wf->Nf + g; 581 582 PetscFunctionBegin; 583 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, i0, (void (*)(void))g0)); 584 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, i1, (void (*)(void))g1)); 585 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, i2, (void (*)(void))g2)); 586 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, i3, (void (*)(void))g3)); 587 PetscFunctionReturn(0); 588 } 589 590 PetscErrorCode PetscWeakFormHasJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre) 591 { 592 PetscInt n0, n1, n2, n3; 593 594 PetscFunctionBegin; 595 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 596 PetscValidBoolPointer(hasJacPre, 2); 597 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP0], &n0)); 598 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP1], &n1)); 599 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP2], &n2)); 600 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP3], &n3)); 601 *hasJacPre = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 602 PetscFunctionReturn(0); 603 } 604 605 PetscErrorCode PetscWeakFormGetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 606 { 607 PetscInt find = f * wf->Nf + g; 608 609 PetscFunctionBegin; 610 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (***)(void))g0)); 611 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (***)(void))g1)); 612 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (***)(void))g2)); 613 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (***)(void))g3)); 614 PetscFunctionReturn(0); 615 } 616 617 PetscErrorCode PetscWeakFormAddJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 618 { 619 PetscInt find = f * wf->Nf + g; 620 621 PetscFunctionBegin; 622 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, (void (*)(void))g0)); 623 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, (void (*)(void))g1)); 624 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, (void (*)(void))g2)); 625 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, (void (*)(void))g3)); 626 PetscFunctionReturn(0); 627 } 628 629 PetscErrorCode PetscWeakFormSetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 630 { 631 PetscInt find = f * wf->Nf + g; 632 633 PetscFunctionBegin; 634 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (**)(void))g0)); 635 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (**)(void))g1)); 636 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (**)(void))g2)); 637 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (**)(void))g3)); 638 PetscFunctionReturn(0); 639 } 640 641 PetscErrorCode PetscWeakFormSetIndexJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 642 { 643 PetscInt find = f * wf->Nf + g; 644 645 PetscFunctionBegin; 646 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, i0, (void (*)(void))g0)); 647 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, i1, (void (*)(void))g1)); 648 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, i2, (void (*)(void))g2)); 649 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, i3, (void (*)(void))g3)); 650 PetscFunctionReturn(0); 651 } 652 653 PetscErrorCode PetscWeakFormHasBdJacobian(PetscWeakForm wf, PetscBool *hasJac) 654 { 655 PetscInt n0, n1, n2, n3; 656 657 PetscFunctionBegin; 658 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 659 PetscValidBoolPointer(hasJac, 2); 660 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG0], &n0)); 661 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG1], &n1)); 662 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG2], &n2)); 663 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG3], &n3)); 664 *hasJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 665 PetscFunctionReturn(0); 666 } 667 668 PetscErrorCode PetscWeakFormGetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 669 { 670 PetscInt find = f * wf->Nf + g; 671 672 PetscFunctionBegin; 673 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (***)(void))g0)); 674 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (***)(void))g1)); 675 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (***)(void))g2)); 676 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (***)(void))g3)); 677 PetscFunctionReturn(0); 678 } 679 680 PetscErrorCode PetscWeakFormAddBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 681 { 682 PetscInt find = f * wf->Nf + g; 683 684 PetscFunctionBegin; 685 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, (void (*)(void))g0)); 686 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, (void (*)(void))g1)); 687 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, (void (*)(void))g2)); 688 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, (void (*)(void))g3)); 689 PetscFunctionReturn(0); 690 } 691 692 PetscErrorCode PetscWeakFormSetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 693 { 694 PetscInt find = f * wf->Nf + g; 695 696 PetscFunctionBegin; 697 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (**)(void))g0)); 698 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (**)(void))g1)); 699 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (**)(void))g2)); 700 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (**)(void))g3)); 701 PetscFunctionReturn(0); 702 } 703 704 PetscErrorCode PetscWeakFormSetIndexBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 705 { 706 PetscInt find = f * wf->Nf + g; 707 708 PetscFunctionBegin; 709 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, i0, (void (*)(void))g0)); 710 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, i1, (void (*)(void))g1)); 711 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, i2, (void (*)(void))g2)); 712 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, i3, (void (*)(void))g3)); 713 PetscFunctionReturn(0); 714 } 715 716 PetscErrorCode PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre) 717 { 718 PetscInt n0, n1, n2, n3; 719 720 PetscFunctionBegin; 721 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 722 PetscValidBoolPointer(hasJacPre, 2); 723 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP0], &n0)); 724 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP1], &n1)); 725 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP2], &n2)); 726 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP3], &n3)); 727 *hasJacPre = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 728 PetscFunctionReturn(0); 729 } 730 731 PetscErrorCode PetscWeakFormGetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 732 { 733 PetscInt find = f * wf->Nf + g; 734 735 PetscFunctionBegin; 736 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (***)(void))g0)); 737 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (***)(void))g1)); 738 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (***)(void))g2)); 739 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (***)(void))g3)); 740 PetscFunctionReturn(0); 741 } 742 743 PetscErrorCode PetscWeakFormAddBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 744 { 745 PetscInt find = f * wf->Nf + g; 746 747 PetscFunctionBegin; 748 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, (void (*)(void))g0)); 749 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, (void (*)(void))g1)); 750 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, (void (*)(void))g2)); 751 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, (void (*)(void))g3)); 752 PetscFunctionReturn(0); 753 } 754 755 PetscErrorCode PetscWeakFormSetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 756 { 757 PetscInt find = f * wf->Nf + g; 758 759 PetscFunctionBegin; 760 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (**)(void))g0)); 761 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (**)(void))g1)); 762 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (**)(void))g2)); 763 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (**)(void))g3)); 764 PetscFunctionReturn(0); 765 } 766 767 PetscErrorCode PetscWeakFormSetIndexBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 768 { 769 PetscInt find = f * wf->Nf + g; 770 771 PetscFunctionBegin; 772 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, i0, (void (*)(void))g0)); 773 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, i1, (void (*)(void))g1)); 774 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, i2, (void (*)(void))g2)); 775 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, i3, (void (*)(void))g3)); 776 PetscFunctionReturn(0); 777 } 778 779 PetscErrorCode PetscWeakFormHasDynamicJacobian(PetscWeakForm wf, PetscBool *hasDynJac) 780 { 781 PetscInt n0, n1, n2, n3; 782 783 PetscFunctionBegin; 784 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 785 PetscValidBoolPointer(hasDynJac, 2); 786 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT0], &n0)); 787 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT1], &n1)); 788 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT2], &n2)); 789 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT3], &n3)); 790 *hasDynJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 791 PetscFunctionReturn(0); 792 } 793 794 PetscErrorCode PetscWeakFormGetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 795 { 796 PetscInt find = f * wf->Nf + g; 797 798 PetscFunctionBegin; 799 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (***)(void))g0)); 800 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (***)(void))g1)); 801 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (***)(void))g2)); 802 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (***)(void))g3)); 803 PetscFunctionReturn(0); 804 } 805 806 PetscErrorCode PetscWeakFormAddDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 807 { 808 PetscInt find = f * wf->Nf + g; 809 810 PetscFunctionBegin; 811 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, (void (*)(void))g0)); 812 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, (void (*)(void))g1)); 813 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, (void (*)(void))g2)); 814 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, (void (*)(void))g3)); 815 PetscFunctionReturn(0); 816 } 817 818 PetscErrorCode PetscWeakFormSetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 819 { 820 PetscInt find = f * wf->Nf + g; 821 822 PetscFunctionBegin; 823 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (**)(void))g0)); 824 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (**)(void))g1)); 825 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (**)(void))g2)); 826 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (**)(void))g3)); 827 PetscFunctionReturn(0); 828 } 829 830 PetscErrorCode PetscWeakFormSetIndexDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 831 { 832 PetscInt find = f * wf->Nf + g; 833 834 PetscFunctionBegin; 835 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, i0, (void (*)(void))g0)); 836 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, i1, (void (*)(void))g1)); 837 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, i2, (void (*)(void))g2)); 838 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, i3, (void (*)(void))g3)); 839 PetscFunctionReturn(0); 840 } 841 842 PetscErrorCode PetscWeakFormGetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n, void (***r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)) 843 { 844 PetscFunctionBegin; 845 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (***)(void))r)); 846 PetscFunctionReturn(0); 847 } 848 849 PetscErrorCode PetscWeakFormSetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n, void (**r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)) 850 { 851 PetscFunctionBegin; 852 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (**)(void))r)); 853 PetscFunctionReturn(0); 854 } 855 856 PetscErrorCode PetscWeakFormSetIndexRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt i, void (*r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)) 857 { 858 PetscFunctionBegin; 859 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, i, (void (*)(void))r)); 860 PetscFunctionReturn(0); 861 } 862 863 /*@ 864 PetscWeakFormGetNumFields - Returns the number of fields in a `PetscWeakForm` 865 866 Not collective 867 868 Input Parameter: 869 . wf - The `PetscWeakForm` object 870 871 Output Parameter: 872 . Nf - The number of fields 873 874 Level: beginner 875 876 .seealso: `PetscWeakForm`, `PetscWeakFormSetNumFields()`, `PetscWeakFormCreate()` 877 @*/ 878 PetscErrorCode PetscWeakFormGetNumFields(PetscWeakForm wf, PetscInt *Nf) 879 { 880 PetscFunctionBegin; 881 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 882 PetscValidIntPointer(Nf, 2); 883 *Nf = wf->Nf; 884 PetscFunctionReturn(0); 885 } 886 887 /*@ 888 PetscWeakFormSetNumFields - Sets the number of fields 889 890 Not collective 891 892 Input Parameters: 893 + wf - The `PetscWeakForm` object 894 - Nf - The number of fields 895 896 Level: beginner 897 898 .seealso: `PetscWeakForm`, `PetscWeakFormGetNumFields()`, `PetscWeakFormCreate()` 899 @*/ 900 PetscErrorCode PetscWeakFormSetNumFields(PetscWeakForm wf, PetscInt Nf) 901 { 902 PetscFunctionBegin; 903 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 904 wf->Nf = Nf; 905 PetscFunctionReturn(0); 906 } 907 908 /*@ 909 PetscWeakFormDestroy - Destroys a `PetscWeakForm` object 910 911 Collective on wf 912 913 Input Parameter: 914 . wf - the `PetscWeakForm` object to destroy 915 916 Level: developer 917 918 .seealso: `PetscWeakForm`, `PetscWeakFormCreate()`, `PetscWeakFormView()` 919 @*/ 920 PetscErrorCode PetscWeakFormDestroy(PetscWeakForm *wf) 921 { 922 PetscInt f; 923 924 PetscFunctionBegin; 925 if (!*wf) PetscFunctionReturn(0); 926 PetscValidHeaderSpecific((*wf), PETSCWEAKFORM_CLASSID, 1); 927 928 if (--((PetscObject)(*wf))->refct > 0) { 929 *wf = NULL; 930 PetscFunctionReturn(0); 931 } 932 ((PetscObject)(*wf))->refct = 0; 933 PetscCall(PetscChunkBufferDestroy(&(*wf)->funcs)); 934 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormDestroy(&(*wf)->form[f])); 935 PetscCall(PetscFree((*wf)->form)); 936 PetscCall(PetscHeaderDestroy(wf)); 937 PetscFunctionReturn(0); 938 } 939 940 static PetscErrorCode PetscWeakFormViewTable_Ascii(PetscWeakForm wf, PetscViewer viewer, PetscBool splitField, const char tableName[], PetscHMapForm map) 941 { 942 PetscInt Nf = wf->Nf, Nk, k; 943 944 PetscFunctionBegin; 945 PetscCall(PetscHMapFormGetSize(map, &Nk)); 946 if (Nk) { 947 PetscFormKey *keys; 948 void (**funcs)(void) = NULL; 949 const char **names; 950 PetscInt *values, *idx1, *idx2, *idx; 951 PetscBool showPart = PETSC_FALSE, showPointer = PETSC_FALSE; 952 PetscInt off = 0; 953 954 PetscCall(PetscMalloc6(Nk, &keys, Nk, &names, Nk, &values, Nk, &idx1, Nk, &idx2, Nk, &idx)); 955 PetscCall(PetscHMapFormGetKeys(map, &off, keys)); 956 /* Sort keys by label name and value */ 957 { 958 /* First sort values */ 959 for (k = 0; k < Nk; ++k) { 960 values[k] = keys[k].value; 961 idx1[k] = k; 962 } 963 PetscCall(PetscSortIntWithPermutation(Nk, values, idx1)); 964 /* If the string sort is stable, it will be sorted correctly overall */ 965 for (k = 0; k < Nk; ++k) { 966 if (keys[idx1[k]].label) PetscCall(PetscObjectGetName((PetscObject)keys[idx1[k]].label, &names[k])); 967 else names[k] = ""; 968 idx2[k] = k; 969 } 970 PetscCall(PetscSortStrWithPermutation(Nk, names, idx2)); 971 for (k = 0; k < Nk; ++k) { 972 if (keys[k].label) PetscCall(PetscObjectGetName((PetscObject)keys[k].label, &names[k])); 973 else names[k] = ""; 974 idx[k] = idx1[idx2[k]]; 975 } 976 } 977 PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", tableName)); 978 PetscCall(PetscViewerASCIIPushTab(viewer)); 979 for (k = 0; k < Nk; ++k) { 980 if (keys[k].part != 0) showPart = PETSC_TRUE; 981 } 982 for (k = 0; k < Nk; ++k) { 983 const PetscInt i = idx[k]; 984 PetscInt n, f; 985 986 if (keys[i].label) { 987 if (showPointer) PetscCall(PetscViewerASCIIPrintf(viewer, "(%s:%p, %" PetscInt_FMT ") ", names[i], (void *)keys[i].label, keys[i].value)); 988 else PetscCall(PetscViewerASCIIPrintf(viewer, "(%s, %" PetscInt_FMT ") ", names[i], keys[i].value)); 989 } 990 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 991 if (splitField) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ", %" PetscInt_FMT ") ", keys[i].field / Nf, keys[i].field % Nf)); 992 else PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") ", keys[i].field)); 993 if (showPart) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") ", keys[i].part)); 994 PetscCall(PetscWeakFormGetFunction_Private(wf, map, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &n, &funcs)); 995 for (f = 0; f < n; ++f) { 996 char *fname; 997 size_t len, l; 998 999 if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 1000 PetscCall(PetscDLAddr(funcs[f], &fname)); 1001 if (fname) { 1002 /* Eliminate argument types */ 1003 PetscCall(PetscStrlen(fname, &len)); 1004 for (l = 0; l < len; ++l) 1005 if (fname[l] == '(') { 1006 fname[l] = '\0'; 1007 break; 1008 } 1009 PetscCall(PetscViewerASCIIPrintf(viewer, "%s", fname)); 1010 } else if (showPointer) { 1011 #if defined(__clang__) 1012 #pragma clang diagnostic push 1013 #pragma clang diagnostic ignored "-Wformat-pedantic" 1014 #elif defined(__GNUC__) || defined(__GNUG__) 1015 #pragma GCC diagnostic push 1016 #pragma GCC diagnostic ignored "-Wformat" 1017 #endif 1018 PetscCall(PetscViewerASCIIPrintf(viewer, "%p", funcs[f])); 1019 #if defined(__clang__) 1020 #pragma clang diagnostic pop 1021 #elif defined(__GNUC__) || defined(__GNUG__) 1022 #pragma GCC diagnostic pop 1023 #endif 1024 } 1025 PetscCall(PetscFree(fname)); 1026 } 1027 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 1028 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 1029 } 1030 PetscCall(PetscViewerASCIIPopTab(viewer)); 1031 PetscCall(PetscFree6(keys, names, values, idx1, idx2, idx)); 1032 } 1033 PetscFunctionReturn(0); 1034 } 1035 1036 static PetscErrorCode PetscWeakFormView_Ascii(PetscWeakForm wf, PetscViewer viewer) 1037 { 1038 PetscViewerFormat format; 1039 PetscInt f; 1040 1041 PetscFunctionBegin; 1042 PetscCall(PetscViewerGetFormat(viewer, &format)); 1043 PetscCall(PetscViewerASCIIPrintf(viewer, "Weak Form System with %" PetscInt_FMT " fields\n", wf->Nf)); 1044 PetscCall(PetscViewerASCIIPushTab(viewer)); 1045 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormViewTable_Ascii(wf, viewer, PETSC_TRUE, PetscWeakFormKinds[f], wf->form[f])); 1046 PetscCall(PetscViewerASCIIPopTab(viewer)); 1047 PetscFunctionReturn(0); 1048 } 1049 1050 /*@C 1051 PetscWeakFormView - Views a `PetscWeakForm` 1052 1053 Collective on wf 1054 1055 Input Parameters: 1056 + wf - the `PetscWeakForm` object to view 1057 - v - the viewer 1058 1059 Level: developer 1060 1061 .seealso: `PetscViewer`, `PetscWeakForm`, `PetscWeakFormDestroy()`, `PetscWeakFormCreate()` 1062 @*/ 1063 PetscErrorCode PetscWeakFormView(PetscWeakForm wf, PetscViewer v) 1064 { 1065 PetscBool iascii; 1066 1067 PetscFunctionBegin; 1068 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 1069 if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)wf), &v)); 1070 else PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2); 1071 PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii)); 1072 if (iascii) PetscCall(PetscWeakFormView_Ascii(wf, v)); 1073 PetscTryTypeMethod(wf, view, v); 1074 PetscFunctionReturn(0); 1075 } 1076 1077 /*@ 1078 PetscWeakFormCreate - Creates an empty `PetscWeakForm` object. 1079 1080 Collective 1081 1082 Input Parameter: 1083 . comm - The communicator for the `PetscWeakForm` object 1084 1085 Output Parameter: 1086 . wf - The `PetscWeakForm` object 1087 1088 Level: beginner 1089 1090 .seealso: `PetscWeakForm`, `PetscDS`, `PetscWeakFormDestroy()` 1091 @*/ 1092 PetscErrorCode PetscWeakFormCreate(MPI_Comm comm, PetscWeakForm *wf) 1093 { 1094 PetscWeakForm p; 1095 PetscInt f; 1096 1097 PetscFunctionBegin; 1098 PetscValidPointer(wf, 2); 1099 *wf = NULL; 1100 PetscCall(PetscDSInitializePackage()); 1101 1102 PetscCall(PetscHeaderCreate(p, PETSCWEAKFORM_CLASSID, "PetscWeakForm", "Weak Form System", "PetscWeakForm", comm, PetscWeakFormDestroy, PetscWeakFormView)); 1103 1104 p->Nf = 0; 1105 PetscCall(PetscChunkBufferCreate(sizeof(&PetscWeakFormCreate), 2, &p->funcs)); 1106 PetscCall(PetscMalloc1(PETSC_NUM_WF, &p->form)); 1107 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormCreate(&p->form[f])); 1108 *wf = p; 1109 PetscFunctionReturn(0); 1110 } 1111