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