16528b96dSMatthew G. Knepley #include <petsc/private/petscdsimpl.h> /*I "petscds.h" I*/ 26528b96dSMatthew G. Knepley 36528b96dSMatthew G. Knepley PetscClassId PETSCWEAKFORM_CLASSID = 0; 46528b96dSMatthew G. Knepley 506ad1575SMatthew G. Knepley 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}; 606ad1575SMatthew G. Knepley 76528b96dSMatthew G. Knepley static PetscErrorCode PetscChunkBufferCreate(size_t unitbytes, size_t expected, PetscChunkBuffer **buffer) 86528b96dSMatthew G. Knepley { 96528b96dSMatthew G. Knepley PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(PetscNew(buffer)); 119566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(expected*unitbytes, &(*buffer)->array)); 126528b96dSMatthew G. Knepley (*buffer)->size = expected; 136528b96dSMatthew G. Knepley (*buffer)->unitbytes = unitbytes; 146528b96dSMatthew G. Knepley (*buffer)->alloc = expected*unitbytes; 156528b96dSMatthew G. Knepley PetscFunctionReturn(0); 166528b96dSMatthew G. Knepley } 176528b96dSMatthew G. Knepley 1845480ffeSMatthew G. Knepley static PetscErrorCode PetscChunkBufferDuplicate(PetscChunkBuffer *buffer, PetscChunkBuffer **bufferNew) 1945480ffeSMatthew G. Knepley { 2045480ffeSMatthew G. Knepley PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(PetscNew(bufferNew)); 229566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(buffer->size*buffer->unitbytes, &(*bufferNew)->array)); 239566063dSJacob Faibussowitsch PetscCall(PetscMemcpy((*bufferNew)->array, buffer->array, buffer->size*buffer->unitbytes)); 2445480ffeSMatthew G. Knepley (*bufferNew)->size = buffer->size; 2545480ffeSMatthew G. Knepley (*bufferNew)->unitbytes = buffer->unitbytes; 2645480ffeSMatthew G. Knepley (*bufferNew)->alloc = buffer->size*buffer->unitbytes; 2745480ffeSMatthew G. Knepley PetscFunctionReturn(0); 2845480ffeSMatthew G. Knepley } 2945480ffeSMatthew G. Knepley 306528b96dSMatthew G. Knepley static PetscErrorCode PetscChunkBufferDestroy(PetscChunkBuffer **buffer) 316528b96dSMatthew G. Knepley { 326528b96dSMatthew G. Knepley PetscFunctionBegin; 339566063dSJacob Faibussowitsch PetscCall(PetscFree((*buffer)->array)); 349566063dSJacob Faibussowitsch PetscCall(PetscFree(*buffer)); 356528b96dSMatthew G. Knepley PetscFunctionReturn(0); 366528b96dSMatthew G. Knepley } 376528b96dSMatthew G. Knepley 386528b96dSMatthew G. Knepley static PetscErrorCode PetscChunkBufferCreateChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk) 396528b96dSMatthew G. Knepley { 406528b96dSMatthew G. Knepley PetscFunctionBegin; 416528b96dSMatthew G. Knepley if ((buffer->size + size)*buffer->unitbytes > buffer->alloc) { 426528b96dSMatthew G. Knepley char *tmp; 436528b96dSMatthew G. Knepley 446528b96dSMatthew G. Knepley if (!buffer->alloc) buffer->alloc = (buffer->size + size)*buffer->unitbytes; 456528b96dSMatthew G. Knepley while ((buffer->size + size)*buffer->unitbytes > buffer->alloc) buffer->alloc *= 2; 469566063dSJacob Faibussowitsch PetscCall(PetscMalloc(buffer->alloc, &tmp)); 479566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(tmp, buffer->array, buffer->size*buffer->unitbytes)); 489566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer->array)); 496528b96dSMatthew G. Knepley buffer->array = tmp; 506528b96dSMatthew G. Knepley } 512baeeceeSMatthew G. Knepley chunk->start = buffer->size*buffer->unitbytes; 526528b96dSMatthew G. Knepley chunk->size = size; 536528b96dSMatthew G. Knepley chunk->reserved = size; 546528b96dSMatthew G. Knepley buffer->size += size; 556528b96dSMatthew G. Knepley PetscFunctionReturn(0); 566528b96dSMatthew G. Knepley } 576528b96dSMatthew G. Knepley 586528b96dSMatthew G. Knepley static PetscErrorCode PetscChunkBufferEnlargeChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk) 596528b96dSMatthew G. Knepley { 606528b96dSMatthew G. Knepley size_t siz = size; 616528b96dSMatthew G. Knepley 626528b96dSMatthew G. Knepley PetscFunctionBegin; 636528b96dSMatthew G. Knepley if (chunk->size + size > chunk->reserved) { 646528b96dSMatthew G. Knepley PetscChunk newchunk; 656528b96dSMatthew G. Knepley PetscInt reserved = chunk->size; 666528b96dSMatthew G. Knepley 676528b96dSMatthew G. Knepley /* TODO Here if we had a chunk list, we could update them all to reclaim unused space */ 686528b96dSMatthew G. Knepley while (reserved < chunk->size+size) reserved *= 2; 699566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(buffer, (size_t) reserved, &newchunk)); 706528b96dSMatthew G. Knepley newchunk.size = chunk->size+size; 719566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(&buffer->array[newchunk.start], &buffer->array[chunk->start], chunk->size * buffer->unitbytes)); 726528b96dSMatthew G. Knepley *chunk = newchunk; 736528b96dSMatthew G. Knepley } else { 746528b96dSMatthew G. Knepley chunk->size += siz; 756528b96dSMatthew G. Knepley } 766528b96dSMatthew G. Knepley PetscFunctionReturn(0); 776528b96dSMatthew G. Knepley } 786528b96dSMatthew G. Knepley 796528b96dSMatthew G. Knepley /*@C 8006ad1575SMatthew G. Knepley PetscFormKeySort - Sorts an array of PetscFormKey in place in increasing order. 816528b96dSMatthew G. Knepley 826528b96dSMatthew G. Knepley Not Collective 836528b96dSMatthew G. Knepley 846528b96dSMatthew G. Knepley Input Parameters: 856528b96dSMatthew G. Knepley + n - number of values 8606ad1575SMatthew G. Knepley - X - array of PetscFormKey 876528b96dSMatthew G. Knepley 886528b96dSMatthew G. Knepley Level: intermediate 896528b96dSMatthew G. Knepley 90db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()` 916528b96dSMatthew G. Knepley @*/ 9206ad1575SMatthew G. Knepley PetscErrorCode PetscFormKeySort(PetscInt n, PetscFormKey arr[]) 936528b96dSMatthew G. Knepley { 946528b96dSMatthew G. Knepley PetscFunctionBegin; 956528b96dSMatthew G. Knepley if (n <= 1) PetscFunctionReturn(0); 966528b96dSMatthew G. Knepley PetscValidPointer(arr, 2); 979566063dSJacob Faibussowitsch PetscCall(PetscTimSort(n, arr, sizeof(PetscFormKey), Compare_PetscFormKey_Private, NULL)); 986528b96dSMatthew G. Knepley PetscFunctionReturn(0); 996528b96dSMatthew G. Knepley } 1006528b96dSMatthew G. Knepley 10106ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt *n, void (***func)()) 1026528b96dSMatthew G. Knepley { 10306ad1575SMatthew G. Knepley PetscFormKey key; 1046528b96dSMatthew G. Knepley PetscChunk chunk; 1056528b96dSMatthew G. Knepley 1066528b96dSMatthew G. Knepley PetscFunctionBegin; 10706ad1575SMatthew G. Knepley key.label = label; key.value = value; key.field = f; key.part = part; 1089566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1096528b96dSMatthew G. Knepley if (chunk.size < 0) {*n = 0; *func = NULL;} 1102baeeceeSMatthew G. Knepley else {*n = chunk.size; *func = (void (**)()) &wf->funcs->array[chunk.start];} 1116528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1126528b96dSMatthew G. Knepley } 1136528b96dSMatthew G. Knepley 1146528b96dSMatthew G. Knepley /* A NULL argument for func causes this to clear the key */ 11506ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt n, void (**func)()) 1166528b96dSMatthew G. Knepley { 11706ad1575SMatthew G. Knepley PetscFormKey key; 1186528b96dSMatthew G. Knepley PetscChunk chunk; 1196528b96dSMatthew G. Knepley PetscInt i; 1206528b96dSMatthew G. Knepley 1216528b96dSMatthew G. Knepley PetscFunctionBegin; 12206ad1575SMatthew G. Knepley key.label = label; key.value = value; key.field = f; key.part = part; 1236528b96dSMatthew G. Knepley if (!func) { 1249566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDel(ht, key)); 1256528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1261baa6e33SBarry Smith } else PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1276528b96dSMatthew G. Knepley if (chunk.size < 0) { 1289566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(wf->funcs, n, &chunk)); 1299566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1306528b96dSMatthew G. Knepley } else if (chunk.size <= n) { 1319566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, n - chunk.size, &chunk)); 1329566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1336528b96dSMatthew G. Knepley } 1342baeeceeSMatthew G. Knepley for (i = 0; i < n; ++i) ((void (**)()) &wf->funcs->array[chunk.start])[i] = func[i]; 1356528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1366528b96dSMatthew G. Knepley } 1376528b96dSMatthew G. Knepley 13806ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormAddFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, void (*func)()) 1396528b96dSMatthew G. Knepley { 14006ad1575SMatthew G. Knepley PetscFormKey key; 1416528b96dSMatthew G. Knepley PetscChunk chunk; 1426528b96dSMatthew G. Knepley 1436528b96dSMatthew G. Knepley PetscFunctionBegin; 1446528b96dSMatthew G. Knepley if (!func) PetscFunctionReturn(0); 14506ad1575SMatthew G. Knepley key.label = label; key.value = value; key.field = f; key.part = part; 1469566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1476528b96dSMatthew G. Knepley if (chunk.size < 0) { 1489566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(wf->funcs, 1, &chunk)); 1499566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1502baeeceeSMatthew G. Knepley ((void (**)()) &wf->funcs->array[chunk.start])[0] = func; 1516528b96dSMatthew G. Knepley } else { 1529566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, 1, &chunk)); 1539566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1542baeeceeSMatthew G. Knepley ((void (**)()) &wf->funcs->array[chunk.start])[chunk.size-1] = func; 1556528b96dSMatthew G. Knepley } 1566528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1576528b96dSMatthew G. Knepley } 1586528b96dSMatthew G. Knepley 15906ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, void (**func)()) 1606528b96dSMatthew G. Knepley { 16106ad1575SMatthew G. Knepley PetscFormKey key; 1626528b96dSMatthew G. Knepley PetscChunk chunk; 1636528b96dSMatthew G. Knepley 1646528b96dSMatthew G. Knepley PetscFunctionBegin; 16506ad1575SMatthew G. Knepley key.label = label; key.value = value; key.field = f; key.part = part; 1669566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1676528b96dSMatthew G. Knepley if (chunk.size < 0) {*func = NULL;} 1686528b96dSMatthew G. Knepley else { 16963a3b9bcSJacob Faibussowitsch PetscCheck(ind < chunk.size,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", ind, chunk.size); 1702baeeceeSMatthew G. Knepley *func = ((void (**)()) &wf->funcs->array[chunk.start])[ind]; 1716528b96dSMatthew G. Knepley } 1726528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1736528b96dSMatthew G. Knepley } 1746528b96dSMatthew G. Knepley 17506ad1575SMatthew G. Knepley /* Ignore a NULL func */ 17606ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, void (*func)()) 1776528b96dSMatthew G. Knepley { 17806ad1575SMatthew G. Knepley PetscFormKey key; 1796528b96dSMatthew G. Knepley PetscChunk chunk; 1806528b96dSMatthew G. Knepley 1816528b96dSMatthew G. Knepley PetscFunctionBegin; 18206ad1575SMatthew G. Knepley if (!func) PetscFunctionReturn(0); 18306ad1575SMatthew G. Knepley key.label = label; key.value = value; key.field = f; key.part = part; 1849566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1856528b96dSMatthew G. Knepley if (chunk.size < 0) { 1869566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(wf->funcs, ind+1, &chunk)); 1879566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1886528b96dSMatthew G. Knepley } else if (chunk.size <= ind) { 1899566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, ind - chunk.size + 1, &chunk)); 1909566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1916528b96dSMatthew G. Knepley } 1922baeeceeSMatthew G. Knepley ((void (**)()) &wf->funcs->array[chunk.start])[ind] = func; 1936528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1946528b96dSMatthew G. Knepley } 1956528b96dSMatthew G. Knepley 19606ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormClearIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind) 19706ad1575SMatthew G. Knepley { 19806ad1575SMatthew G. Knepley PetscFormKey key; 19906ad1575SMatthew G. Knepley PetscChunk chunk; 20006ad1575SMatthew G. Knepley 20106ad1575SMatthew G. Knepley PetscFunctionBegin; 20206ad1575SMatthew G. Knepley key.label = label; key.value = value; key.field = f; key.part = part; 2039566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 20406ad1575SMatthew G. Knepley if (chunk.size < 0) { 20506ad1575SMatthew G. Knepley PetscFunctionReturn(0); 20606ad1575SMatthew G. Knepley } else if (!ind && chunk.size == 1) { 2079566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDel(ht, key)); 20806ad1575SMatthew G. Knepley PetscFunctionReturn(0); 20906ad1575SMatthew G. Knepley } else if (chunk.size <= ind) { 21006ad1575SMatthew G. Knepley PetscFunctionReturn(0); 21106ad1575SMatthew G. Knepley } 21206ad1575SMatthew G. Knepley ((void (**)()) &wf->funcs->array[chunk.start])[ind] = NULL; 21306ad1575SMatthew G. Knepley PetscFunctionReturn(0); 21406ad1575SMatthew G. Knepley } 21506ad1575SMatthew G. Knepley 21645480ffeSMatthew G. Knepley /*@ 21745480ffeSMatthew G. Knepley PetscWeakFormCopy - Copy the pointwise functions to another PetscWeakForm 21845480ffeSMatthew G. Knepley 21945480ffeSMatthew G. Knepley Not Collective 22045480ffeSMatthew G. Knepley 22145480ffeSMatthew G. Knepley Input Parameter: 22245480ffeSMatthew G. Knepley . wf - The original PetscWeakForm 22345480ffeSMatthew G. Knepley 22445480ffeSMatthew G. Knepley Output Parameter: 22545480ffeSMatthew G. Knepley . wfNew - The copy PetscWeakForm 22645480ffeSMatthew G. Knepley 22745480ffeSMatthew G. Knepley Level: intermediate 22845480ffeSMatthew G. Knepley 229db781477SPatrick Sanan .seealso: `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 23045480ffeSMatthew G. Knepley @*/ 23145480ffeSMatthew G. Knepley PetscErrorCode PetscWeakFormCopy(PetscWeakForm wf, PetscWeakForm wfNew) 23245480ffeSMatthew G. Knepley { 23306ad1575SMatthew G. Knepley PetscInt f; 23445480ffeSMatthew G. Knepley 23545480ffeSMatthew G. Knepley PetscFunctionBegin; 23645480ffeSMatthew G. Knepley wfNew->Nf = wf->Nf; 2379566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferDestroy(&wfNew->funcs)); 2389566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferDuplicate(wf->funcs, &wfNew->funcs)); 23906ad1575SMatthew G. Knepley for (f = 0; f < PETSC_NUM_WF; ++f) { 2409566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDestroy(&wfNew->form[f])); 2419566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDuplicate(wf->form[f], &wfNew->form[f])); 24206ad1575SMatthew G. Knepley } 24306ad1575SMatthew G. Knepley PetscFunctionReturn(0); 24406ad1575SMatthew G. Knepley } 24506ad1575SMatthew G. Knepley 24606ad1575SMatthew G. Knepley /*@ 24706ad1575SMatthew G. Knepley PetscWeakFormClear - Clear all functions from the PetscWeakForm 24806ad1575SMatthew G. Knepley 24906ad1575SMatthew G. Knepley Not Collective 25006ad1575SMatthew G. Knepley 25106ad1575SMatthew G. Knepley Input Parameter: 25206ad1575SMatthew G. Knepley . wf - The original PetscWeakForm 25306ad1575SMatthew G. Knepley 25406ad1575SMatthew G. Knepley Level: intermediate 25506ad1575SMatthew G. Knepley 256db781477SPatrick Sanan .seealso: `PetscWeakFormCopy()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 25706ad1575SMatthew G. Knepley @*/ 25806ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormClear(PetscWeakForm wf) 25906ad1575SMatthew G. Knepley { 26006ad1575SMatthew G. Knepley PetscInt f; 26106ad1575SMatthew G. Knepley 26206ad1575SMatthew G. Knepley PetscFunctionBegin; 2639566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormClear(wf->form[f])); 26445480ffeSMatthew G. Knepley PetscFunctionReturn(0); 26545480ffeSMatthew G. Knepley } 26645480ffeSMatthew G. Knepley 26745480ffeSMatthew G. Knepley static PetscErrorCode PetscWeakFormRewriteKeys_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label, PetscInt Nv, const PetscInt values[]) 26845480ffeSMatthew G. Knepley { 26906ad1575SMatthew G. Knepley PetscFormKey *keys; 27045480ffeSMatthew G. Knepley PetscInt n, i, v, off = 0; 27145480ffeSMatthew G. Knepley 27245480ffeSMatthew G. Knepley PetscFunctionBegin; 2739566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(hmap, &n)); 2749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &keys)); 2759566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetKeys(hmap, &off, keys)); 27645480ffeSMatthew G. Knepley for (i = 0; i < n; ++i) { 27745480ffeSMatthew G. Knepley if (keys[i].label == label) { 27845480ffeSMatthew G. Knepley PetscBool clear = PETSC_TRUE; 27945480ffeSMatthew G. Knepley void (**funcs)(); 28045480ffeSMatthew G. Knepley PetscInt Nf; 28145480ffeSMatthew G. Knepley 2829566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 28345480ffeSMatthew G. Knepley for (v = 0; v < Nv; ++v) { 2849566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, values[v], keys[i].field, keys[i].part, Nf, funcs)); 28545480ffeSMatthew G. Knepley if (values[v] == keys[i].value) clear = PETSC_FALSE; 28645480ffeSMatthew G. Knepley } 2879566063dSJacob Faibussowitsch if (clear) PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL)); 28845480ffeSMatthew G. Knepley } 28945480ffeSMatthew G. Knepley } 2909566063dSJacob Faibussowitsch PetscCall(PetscFree(keys)); 29145480ffeSMatthew G. Knepley PetscFunctionReturn(0); 29245480ffeSMatthew G. Knepley } 29345480ffeSMatthew G. Knepley 29445480ffeSMatthew G. Knepley /*@C 29545480ffeSMatthew G. Knepley PetscWeakFormRewriteKeys - Change any key on the given label to use the new set of label values 29645480ffeSMatthew G. Knepley 29745480ffeSMatthew G. Knepley Not Collective 29845480ffeSMatthew G. Knepley 29945480ffeSMatthew G. Knepley Input Parameters: 30045480ffeSMatthew G. Knepley + wf - The original PetscWeakForm 30145480ffeSMatthew G. Knepley . label - The label to change keys for 30245480ffeSMatthew G. Knepley . Nv - The number of new label values 30345480ffeSMatthew G. Knepley - values - The set of new values to relabel keys with 30445480ffeSMatthew G. Knepley 30545480ffeSMatthew G. Knepley Note: This is used internally when boundary label values are specified from the command line. 30645480ffeSMatthew G. Knepley 30745480ffeSMatthew G. Knepley Level: intermediate 30845480ffeSMatthew G. Knepley 309db781477SPatrick Sanan .seealso: `PetscWeakFormReplaceLabel()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 31045480ffeSMatthew G. Knepley @*/ 31145480ffeSMatthew G. Knepley PetscErrorCode PetscWeakFormRewriteKeys(PetscWeakForm wf, DMLabel label, PetscInt Nv, const PetscInt values[]) 31245480ffeSMatthew G. Knepley { 31306ad1575SMatthew G. Knepley PetscInt f; 31445480ffeSMatthew G. Knepley 31545480ffeSMatthew G. Knepley PetscFunctionBegin; 3169566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormRewriteKeys_Internal(wf, wf->form[f], label, Nv, values)); 31745480ffeSMatthew G. Knepley PetscFunctionReturn(0); 31845480ffeSMatthew G. Knepley } 31945480ffeSMatthew G. Knepley 320d700741fSMatthew G. Knepley static PetscErrorCode PetscWeakFormReplaceLabel_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label) 321d700741fSMatthew G. Knepley { 322d700741fSMatthew G. Knepley PetscFormKey *keys; 323d700741fSMatthew G. Knepley PetscInt n, i, off = 0, maxFuncs = 0; 324d700741fSMatthew G. Knepley void (**tmpf)(); 325d700741fSMatthew G. Knepley const char *name = NULL; 326d700741fSMatthew G. Knepley 327d700741fSMatthew G. Knepley PetscFunctionBegin; 3289566063dSJacob Faibussowitsch if (label) PetscCall(PetscObjectGetName((PetscObject) label, &name)); 3299566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(hmap, &n)); 3309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &keys)); 3319566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetKeys(hmap, &off, keys)); 332d700741fSMatthew G. Knepley for (i = 0; i < n; ++i) { 333d700741fSMatthew G. Knepley PetscBool match = PETSC_FALSE; 334d700741fSMatthew G. Knepley const char *lname = NULL; 335d700741fSMatthew G. Knepley 336d700741fSMatthew G. Knepley if (label == keys[i].label) continue; 3379566063dSJacob Faibussowitsch if (keys[i].label) PetscCall(PetscObjectGetName((PetscObject) keys[i].label, &lname)); 3389566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(name, lname, &match)); 339d700741fSMatthew G. Knepley if ((!name && !lname) || match) { 340d700741fSMatthew G. Knepley void (**funcs)(); 341d700741fSMatthew G. Knepley PetscInt Nf; 342d700741fSMatthew G. Knepley 3439566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 344d700741fSMatthew G. Knepley maxFuncs = PetscMax(maxFuncs, Nf); 345d700741fSMatthew G. Knepley } 346d700741fSMatthew G. Knepley } 347d700741fSMatthew G. Knepley /* Need temp space because chunk buffer can be reallocated in SetFunction() call */ 3489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxFuncs, &tmpf)); 349d700741fSMatthew G. Knepley for (i = 0; i < n; ++i) { 350d700741fSMatthew G. Knepley PetscBool match = PETSC_FALSE; 351d700741fSMatthew G. Knepley const char *lname = NULL; 352d700741fSMatthew G. Knepley 353d700741fSMatthew G. Knepley if (label == keys[i].label) continue; 3549566063dSJacob Faibussowitsch if (keys[i].label) PetscCall(PetscObjectGetName((PetscObject) keys[i].label, &lname)); 3559566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(name, lname, &match)); 356d700741fSMatthew G. Knepley if ((!name && !lname) || match) { 357d700741fSMatthew G. Knepley void (**funcs)(); 358d700741fSMatthew G. Knepley PetscInt Nf, j; 359d700741fSMatthew G. Knepley 3609566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 361d700741fSMatthew G. Knepley for (j = 0; j < Nf; ++j) tmpf[j] = funcs[j]; 3629566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, label, keys[i].value, keys[i].field, keys[i].part, Nf, tmpf)); 3639566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL)); 364d700741fSMatthew G. Knepley } 365d700741fSMatthew G. Knepley } 3669566063dSJacob Faibussowitsch PetscCall(PetscFree(tmpf)); 3679566063dSJacob Faibussowitsch PetscCall(PetscFree(keys)); 368d700741fSMatthew G. Knepley PetscFunctionReturn(0); 369d700741fSMatthew G. Knepley } 370d700741fSMatthew G. Knepley 371d700741fSMatthew G. Knepley /*@C 372d700741fSMatthew G. Knepley PetscWeakFormReplaceLabel - Change any key on a label of the same name to use the new label 373d700741fSMatthew G. Knepley 374d700741fSMatthew G. Knepley Not Collective 375d700741fSMatthew G. Knepley 376d700741fSMatthew G. Knepley Input Parameters: 377d700741fSMatthew G. Knepley + wf - The original PetscWeakForm 378d700741fSMatthew G. Knepley - label - The label to change keys for 379d700741fSMatthew G. Knepley 380d700741fSMatthew G. Knepley Note: This is used internally when meshes are modified 381d700741fSMatthew G. Knepley 382d700741fSMatthew G. Knepley Level: intermediate 383d700741fSMatthew G. Knepley 384db781477SPatrick Sanan .seealso: `PetscWeakFormRewriteKeys()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 385d700741fSMatthew G. Knepley @*/ 386d700741fSMatthew G. Knepley PetscErrorCode PetscWeakFormReplaceLabel(PetscWeakForm wf, DMLabel label) 387d700741fSMatthew G. Knepley { 388d700741fSMatthew G. Knepley PetscInt f; 389d700741fSMatthew G. Knepley 390d700741fSMatthew G. Knepley PetscFunctionBegin; 3919566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormReplaceLabel_Internal(wf, wf->form[f], label)); 392d700741fSMatthew G. Knepley PetscFunctionReturn(0); 393d700741fSMatthew G. Knepley } 394d700741fSMatthew G. Knepley 39506ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormClearIndex(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscWeakFormKind kind, PetscInt ind) 39606ad1575SMatthew G. Knepley { 39706ad1575SMatthew G. Knepley PetscFunctionBegin; 3989566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClearIndexFunction_Private(wf, wf->form[kind], label, val, f, part, ind)); 39906ad1575SMatthew G. Knepley PetscFunctionReturn(0); 40006ad1575SMatthew G. Knepley } 40106ad1575SMatthew G. Knepley 40206ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n, 4036528b96dSMatthew G. Knepley void (***obj)(PetscInt, PetscInt, PetscInt, 4046528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4056528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4066528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 4076528b96dSMatthew G. Knepley { 4086528b96dSMatthew G. Knepley PetscFunctionBegin; 4099566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (***)(void)) obj)); 4106528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4116528b96dSMatthew G. Knepley } 4126528b96dSMatthew G. Knepley 41306ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n, 4146528b96dSMatthew G. Knepley void (**obj)(PetscInt, PetscInt, PetscInt, 4156528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4166528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4176528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 4186528b96dSMatthew G. Knepley { 4196528b96dSMatthew G. Knepley PetscFunctionBegin; 4209566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (**)(void)) obj)); 4216528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4226528b96dSMatthew G. Knepley } 4236528b96dSMatthew G. Knepley 42406ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormAddObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 4256528b96dSMatthew G. Knepley void (*obj)(PetscInt, PetscInt, PetscInt, 4266528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4276528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4286528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 4296528b96dSMatthew G. Knepley { 4306528b96dSMatthew G. Knepley PetscFunctionBegin; 4319566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, (void (*)(void)) obj)); 4326528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4336528b96dSMatthew G. Knepley } 4346528b96dSMatthew G. Knepley 43506ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt ind, 4366528b96dSMatthew G. Knepley void (**obj)(PetscInt, PetscInt, PetscInt, 4376528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4386528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4396528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 4406528b96dSMatthew G. Knepley { 4416528b96dSMatthew G. Knepley PetscFunctionBegin; 4429566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (void (**)(void)) obj)); 4436528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4446528b96dSMatthew G. Knepley } 4456528b96dSMatthew G. Knepley 44606ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt ind, 4476528b96dSMatthew G. Knepley void (*obj)(PetscInt, PetscInt, PetscInt, 4486528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4496528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4506528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 4516528b96dSMatthew G. Knepley { 4526528b96dSMatthew G. Knepley PetscFunctionBegin; 4539566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (void (*)(void)) obj)); 4546528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4556528b96dSMatthew G. Knepley } 4566528b96dSMatthew G. Knepley 45706ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 4586528b96dSMatthew G. Knepley PetscInt *n0, 4596528b96dSMatthew G. Knepley void (***f0)(PetscInt, PetscInt, PetscInt, 4606528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4616528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4626528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 4636528b96dSMatthew G. Knepley PetscInt *n1, 4646528b96dSMatthew G. Knepley void (***f1)(PetscInt, PetscInt, PetscInt, 4656528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4666528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4676528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 4686528b96dSMatthew G. Knepley { 4696528b96dSMatthew G. Knepley PetscFunctionBegin; 4709566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (***)(void)) f0)); 4719566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (***)(void)) f1)); 4726528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4736528b96dSMatthew G. Knepley } 4746528b96dSMatthew G. Knepley 47506ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormAddResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 4766528b96dSMatthew G. Knepley void (*f0)(PetscInt, PetscInt, PetscInt, 4776528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4786528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4796528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 4806528b96dSMatthew G. Knepley void (*f1)(PetscInt, PetscInt, PetscInt, 4816528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4826528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4836528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 4846528b96dSMatthew G. Knepley { 4856528b96dSMatthew G. Knepley PetscFunctionBegin; 4869566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, (void (*)(void)) f0)); 4879566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, (void (*)(void)) f1)); 4886528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4896528b96dSMatthew G. Knepley } 4906528b96dSMatthew G. Knepley 49106ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 4926528b96dSMatthew G. Knepley PetscInt n0, 4936528b96dSMatthew G. Knepley void (**f0)(PetscInt, PetscInt, PetscInt, 4946528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4956528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4966528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 4976528b96dSMatthew G. Knepley PetscInt n1, 4986528b96dSMatthew G. Knepley void (**f1)(PetscInt, PetscInt, PetscInt, 4996528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5006528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5016528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 5026528b96dSMatthew G. Knepley { 5036528b96dSMatthew G. Knepley PetscFunctionBegin; 5049566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (**)(void)) f0)); 5059566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (**)(void)) f1)); 5066528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5076528b96dSMatthew G. Knepley } 5086528b96dSMatthew G. Knepley 50906ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetIndexResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 5106528b96dSMatthew G. Knepley PetscInt i0, 5116528b96dSMatthew G. Knepley void (*f0)(PetscInt, PetscInt, PetscInt, 5126528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5136528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5146528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 5156528b96dSMatthew G. Knepley PetscInt i1, 5166528b96dSMatthew G. Knepley void (*f1)(PetscInt, PetscInt, PetscInt, 5176528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5186528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5196528b96dSMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 5206528b96dSMatthew G. Knepley { 5216528b96dSMatthew G. Knepley PetscFunctionBegin; 5229566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, i0, (void (*)(void)) f0)); 5239566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, i1, (void (*)(void)) f1)); 5246528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5256528b96dSMatthew G. Knepley } 5266528b96dSMatthew G. Knepley 52706ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 5286528b96dSMatthew G. Knepley PetscInt *n0, 5296528b96dSMatthew G. Knepley void (***f0)(PetscInt, PetscInt, PetscInt, 5306528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5316528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5326528b96dSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 5336528b96dSMatthew G. Knepley PetscInt *n1, 5346528b96dSMatthew G. Knepley void (***f1)(PetscInt, PetscInt, PetscInt, 5356528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5366528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5376528b96dSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 5386528b96dSMatthew G. Knepley { 5396528b96dSMatthew G. Knepley PetscFunctionBegin; 5409566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (***)(void)) f0)); 5419566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (***)(void)) f1)); 5426528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5436528b96dSMatthew G. Knepley } 5446528b96dSMatthew G. Knepley 54506ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormAddBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 5466528b96dSMatthew G. Knepley void (*f0)(PetscInt, PetscInt, PetscInt, 5476528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5486528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5496528b96dSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 5506528b96dSMatthew G. Knepley void (*f1)(PetscInt, PetscInt, PetscInt, 5516528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5526528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5536528b96dSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 5546528b96dSMatthew G. Knepley { 5556528b96dSMatthew G. Knepley PetscFunctionBegin; 5569566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, (void (*)(void)) f0)); 5579566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, (void (*)(void)) f1)); 5586528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5596528b96dSMatthew G. Knepley } 5606528b96dSMatthew G. Knepley 56106ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 5626528b96dSMatthew G. Knepley PetscInt n0, 5636528b96dSMatthew G. Knepley void (**f0)(PetscInt, PetscInt, PetscInt, 5646528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5656528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5666528b96dSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 5676528b96dSMatthew G. Knepley PetscInt n1, 5686528b96dSMatthew G. Knepley void (**f1)(PetscInt, PetscInt, PetscInt, 5696528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5706528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5716528b96dSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 5726528b96dSMatthew G. Knepley { 5736528b96dSMatthew G. Knepley PetscFunctionBegin; 5749566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (**)(void)) f0)); 5759566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (**)(void)) f1)); 5766528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5776528b96dSMatthew G. Knepley } 5786528b96dSMatthew G. Knepley 57906ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetIndexBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 5806528b96dSMatthew G. Knepley PetscInt i0, 5816528b96dSMatthew G. Knepley void (*f0)(PetscInt, PetscInt, PetscInt, 5826528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5836528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5846528b96dSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 5856528b96dSMatthew G. Knepley PetscInt i1, 5866528b96dSMatthew G. Knepley void (*f1)(PetscInt, PetscInt, PetscInt, 5876528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5886528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 5896528b96dSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 5906528b96dSMatthew G. Knepley { 5916528b96dSMatthew G. Knepley PetscFunctionBegin; 5929566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, i0, (void (*)(void)) f0)); 5939566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, i1, (void (*)(void)) f1)); 5946528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5956528b96dSMatthew G. Knepley } 5966528b96dSMatthew G. Knepley 5976528b96dSMatthew G. Knepley PetscErrorCode PetscWeakFormHasJacobian(PetscWeakForm wf, PetscBool *hasJac) 5986528b96dSMatthew G. Knepley { 5996528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 6006528b96dSMatthew G. Knepley 6016528b96dSMatthew G. Knepley PetscFunctionBegin; 6026528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 6036528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJac, 2); 6049566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G0], &n0)); 6059566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G1], &n1)); 6069566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G2], &n2)); 6079566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G3], &n3)); 6086528b96dSMatthew G. Knepley *hasJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE; 6096528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6106528b96dSMatthew G. Knepley } 6116528b96dSMatthew G. Knepley 61206ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 6136528b96dSMatthew G. Knepley PetscInt *n0, 6146528b96dSMatthew G. Knepley void (***g0)(PetscInt, PetscInt, PetscInt, 6156528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6166528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6176528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6186528b96dSMatthew G. Knepley PetscInt *n1, 6196528b96dSMatthew G. Knepley void (***g1)(PetscInt, PetscInt, PetscInt, 6206528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6216528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6226528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6236528b96dSMatthew G. Knepley PetscInt *n2, 6246528b96dSMatthew G. Knepley void (***g2)(PetscInt, PetscInt, PetscInt, 6256528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6266528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6276528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6286528b96dSMatthew G. Knepley PetscInt *n3, 6296528b96dSMatthew G. Knepley void (***g3)(PetscInt, PetscInt, PetscInt, 6306528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6316528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6326528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 6336528b96dSMatthew G. Knepley { 6346528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 6356528b96dSMatthew G. Knepley 6366528b96dSMatthew G. Knepley PetscFunctionBegin; 6379566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (***)(void)) g0)); 6389566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (***)(void)) g1)); 6399566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (***)(void)) g2)); 6409566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (***)(void)) g3)); 6416528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6426528b96dSMatthew G. Knepley } 6436528b96dSMatthew G. Knepley 64406ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormAddJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 6456528b96dSMatthew G. Knepley void (*g0)(PetscInt, PetscInt, PetscInt, 6466528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6476528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6486528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6496528b96dSMatthew G. Knepley void (*g1)(PetscInt, PetscInt, PetscInt, 6506528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6516528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6526528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6536528b96dSMatthew G. Knepley void (*g2)(PetscInt, PetscInt, PetscInt, 6546528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6556528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6566528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6576528b96dSMatthew G. Knepley void (*g3)(PetscInt, PetscInt, PetscInt, 6586528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6596528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6606528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 6616528b96dSMatthew G. Knepley { 6626528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 6636528b96dSMatthew G. Knepley 6646528b96dSMatthew G. Knepley PetscFunctionBegin; 6659566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, (void (*)(void)) g0)); 6669566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, (void (*)(void)) g1)); 6679566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, (void (*)(void)) g2)); 6689566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, (void (*)(void)) g3)); 6696528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6706528b96dSMatthew G. Knepley } 6716528b96dSMatthew G. Knepley 67206ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 6736528b96dSMatthew G. Knepley PetscInt n0, 6746528b96dSMatthew G. Knepley void (**g0)(PetscInt, PetscInt, PetscInt, 6756528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6766528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6776528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6786528b96dSMatthew G. Knepley PetscInt n1, 6796528b96dSMatthew G. Knepley void (**g1)(PetscInt, PetscInt, PetscInt, 6806528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6816528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6826528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6836528b96dSMatthew G. Knepley PetscInt n2, 6846528b96dSMatthew G. Knepley void (**g2)(PetscInt, PetscInt, PetscInt, 6856528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6866528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6876528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6886528b96dSMatthew G. Knepley PetscInt n3, 6896528b96dSMatthew G. Knepley void (**g3)(PetscInt, PetscInt, PetscInt, 6906528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6916528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6926528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 6936528b96dSMatthew G. Knepley { 6946528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 6956528b96dSMatthew G. Knepley 6966528b96dSMatthew G. Knepley PetscFunctionBegin; 6979566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (**)(void)) g0)); 6989566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (**)(void)) g1)); 6999566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (**)(void)) g2)); 7009566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (**)(void)) g3)); 7016528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7026528b96dSMatthew G. Knepley } 7036528b96dSMatthew G. Knepley 70406ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetIndexJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 7056528b96dSMatthew G. Knepley PetscInt i0, 7066528b96dSMatthew G. Knepley void (*g0)(PetscInt, PetscInt, PetscInt, 7076528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7086528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7096528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 7106528b96dSMatthew G. Knepley PetscInt i1, 7116528b96dSMatthew G. Knepley void (*g1)(PetscInt, PetscInt, PetscInt, 7126528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7136528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7146528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 7156528b96dSMatthew G. Knepley PetscInt i2, 7166528b96dSMatthew G. Knepley void (*g2)(PetscInt, PetscInt, PetscInt, 7176528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7186528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7196528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 7206528b96dSMatthew G. Knepley PetscInt i3, 7216528b96dSMatthew G. Knepley void (*g3)(PetscInt, PetscInt, PetscInt, 7226528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7236528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7246528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 7256528b96dSMatthew G. Knepley { 7266528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 7276528b96dSMatthew G. Knepley 7286528b96dSMatthew G. Knepley PetscFunctionBegin; 7299566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, i0, (void (*)(void)) g0)); 7309566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, i1, (void (*)(void)) g1)); 7319566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, i2, (void (*)(void)) g2)); 7329566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, i3, (void (*)(void)) g3)); 7336528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7346528b96dSMatthew G. Knepley } 7356528b96dSMatthew G. Knepley 7366528b96dSMatthew G. Knepley PetscErrorCode PetscWeakFormHasJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre) 7376528b96dSMatthew G. Knepley { 7386528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 7396528b96dSMatthew G. Knepley 7406528b96dSMatthew G. Knepley PetscFunctionBegin; 7416528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 7426528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJacPre, 2); 7439566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP0], &n0)); 7449566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP1], &n1)); 7459566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP2], &n2)); 7469566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP3], &n3)); 7476528b96dSMatthew G. Knepley *hasJacPre = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE; 7486528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7496528b96dSMatthew G. Knepley } 7506528b96dSMatthew G. Knepley 75106ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 7526528b96dSMatthew G. Knepley PetscInt *n0, 7536528b96dSMatthew G. Knepley void (***g0)(PetscInt, PetscInt, PetscInt, 7546528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7556528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7566528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 7576528b96dSMatthew G. Knepley PetscInt *n1, 7586528b96dSMatthew G. Knepley void (***g1)(PetscInt, PetscInt, PetscInt, 7596528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7606528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7616528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 7626528b96dSMatthew G. Knepley PetscInt *n2, 7636528b96dSMatthew G. Knepley void (***g2)(PetscInt, PetscInt, PetscInt, 7646528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7656528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7666528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 7676528b96dSMatthew G. Knepley PetscInt *n3, 7686528b96dSMatthew G. Knepley void (***g3)(PetscInt, PetscInt, PetscInt, 7696528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7706528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7716528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 7726528b96dSMatthew G. Knepley { 7736528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 7746528b96dSMatthew G. Knepley 7756528b96dSMatthew G. Knepley PetscFunctionBegin; 7769566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (***)(void)) g0)); 7779566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (***)(void)) g1)); 7789566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (***)(void)) g2)); 7799566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (***)(void)) g3)); 7806528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7816528b96dSMatthew G. Knepley } 7826528b96dSMatthew G. Knepley 78306ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormAddJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 7846528b96dSMatthew G. Knepley void (*g0)(PetscInt, PetscInt, PetscInt, 7856528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7866528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7876528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 7886528b96dSMatthew G. Knepley void (*g1)(PetscInt, PetscInt, PetscInt, 7896528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7906528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7916528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 7926528b96dSMatthew G. Knepley void (*g2)(PetscInt, PetscInt, PetscInt, 7936528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7946528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7956528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 7966528b96dSMatthew G. Knepley void (*g3)(PetscInt, PetscInt, PetscInt, 7976528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7986528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7996528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 8006528b96dSMatthew G. Knepley { 8016528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 8026528b96dSMatthew G. Knepley 8036528b96dSMatthew G. Knepley PetscFunctionBegin; 8049566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, (void (*)(void)) g0)); 8059566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, (void (*)(void)) g1)); 8069566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, (void (*)(void)) g2)); 8079566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, (void (*)(void)) g3)); 8086528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8096528b96dSMatthew G. Knepley } 8106528b96dSMatthew G. Knepley 81106ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 8126528b96dSMatthew G. Knepley PetscInt n0, 8136528b96dSMatthew G. Knepley void (**g0)(PetscInt, PetscInt, PetscInt, 8146528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8156528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8166528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 8176528b96dSMatthew G. Knepley PetscInt n1, 8186528b96dSMatthew G. Knepley void (**g1)(PetscInt, PetscInt, PetscInt, 8196528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8206528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8216528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 8226528b96dSMatthew G. Knepley PetscInt n2, 8236528b96dSMatthew G. Knepley void (**g2)(PetscInt, PetscInt, PetscInt, 8246528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8256528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8266528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 8276528b96dSMatthew G. Knepley PetscInt n3, 8286528b96dSMatthew G. Knepley void (**g3)(PetscInt, PetscInt, PetscInt, 8296528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8306528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8316528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 8326528b96dSMatthew G. Knepley { 8336528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 8346528b96dSMatthew G. Knepley 8356528b96dSMatthew G. Knepley PetscFunctionBegin; 8369566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (**)(void)) g0)); 8379566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (**)(void)) g1)); 8389566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (**)(void)) g2)); 8399566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (**)(void)) g3)); 8406528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8416528b96dSMatthew G. Knepley } 8426528b96dSMatthew G. Knepley 84306ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetIndexJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 8446528b96dSMatthew G. Knepley PetscInt i0, 8456528b96dSMatthew G. Knepley void (*g0)(PetscInt, PetscInt, PetscInt, 8466528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8476528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8486528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 8496528b96dSMatthew G. Knepley PetscInt i1, 8506528b96dSMatthew G. Knepley void (*g1)(PetscInt, PetscInt, PetscInt, 8516528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8526528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8536528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 8546528b96dSMatthew G. Knepley PetscInt i2, 8556528b96dSMatthew G. Knepley void (*g2)(PetscInt, PetscInt, PetscInt, 8566528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8576528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8586528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 8596528b96dSMatthew G. Knepley PetscInt i3, 8606528b96dSMatthew G. Knepley void (*g3)(PetscInt, PetscInt, PetscInt, 8616528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8626528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8636528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 8646528b96dSMatthew G. Knepley { 8656528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 8666528b96dSMatthew G. Knepley 8676528b96dSMatthew G. Knepley PetscFunctionBegin; 8689566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, i0, (void (*)(void)) g0)); 8699566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, i1, (void (*)(void)) g1)); 8709566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, i2, (void (*)(void)) g2)); 8719566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, i3, (void (*)(void)) g3)); 8726528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8736528b96dSMatthew G. Knepley } 8746528b96dSMatthew G. Knepley 8756528b96dSMatthew G. Knepley PetscErrorCode PetscWeakFormHasBdJacobian(PetscWeakForm wf, PetscBool *hasJac) 8766528b96dSMatthew G. Knepley { 8776528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 8786528b96dSMatthew G. Knepley 8796528b96dSMatthew G. Knepley PetscFunctionBegin; 8806528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 8816528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJac, 2); 8829566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG0], &n0)); 8839566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG1], &n1)); 8849566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG2], &n2)); 8859566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG3], &n3)); 8866528b96dSMatthew G. Knepley *hasJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE; 8876528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8886528b96dSMatthew G. Knepley } 8896528b96dSMatthew G. Knepley 89006ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 8916528b96dSMatthew G. Knepley PetscInt *n0, 8926528b96dSMatthew G. Knepley void (***g0)(PetscInt, PetscInt, PetscInt, 8936528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8946528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8956528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 8966528b96dSMatthew G. Knepley PetscInt *n1, 8976528b96dSMatthew G. Knepley void (***g1)(PetscInt, PetscInt, PetscInt, 8986528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8996528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9006528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9016528b96dSMatthew G. Knepley PetscInt *n2, 9026528b96dSMatthew G. Knepley void (***g2)(PetscInt, PetscInt, PetscInt, 9036528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9046528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9056528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9066528b96dSMatthew G. Knepley PetscInt *n3, 9076528b96dSMatthew G. Knepley void (***g3)(PetscInt, PetscInt, PetscInt, 9086528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9096528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9106528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 9116528b96dSMatthew G. Knepley { 9126528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 9136528b96dSMatthew G. Knepley 9146528b96dSMatthew G. Knepley PetscFunctionBegin; 9159566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (***)(void)) g0)); 9169566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (***)(void)) g1)); 9179566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (***)(void)) g2)); 9189566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (***)(void)) g3)); 9196528b96dSMatthew G. Knepley PetscFunctionReturn(0); 9206528b96dSMatthew G. Knepley } 9216528b96dSMatthew G. Knepley 92206ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormAddBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 9236528b96dSMatthew G. Knepley void (*g0)(PetscInt, PetscInt, PetscInt, 9246528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9256528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9266528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9276528b96dSMatthew G. Knepley void (*g1)(PetscInt, PetscInt, PetscInt, 9286528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9296528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9306528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9316528b96dSMatthew G. Knepley void (*g2)(PetscInt, PetscInt, PetscInt, 9326528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9336528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9346528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9356528b96dSMatthew G. Knepley void (*g3)(PetscInt, PetscInt, PetscInt, 9366528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9376528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9386528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 9396528b96dSMatthew G. Knepley { 9406528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 9416528b96dSMatthew G. Knepley 9426528b96dSMatthew G. Knepley PetscFunctionBegin; 9439566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, (void (*)(void)) g0)); 9449566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, (void (*)(void)) g1)); 9459566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, (void (*)(void)) g2)); 9469566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, (void (*)(void)) g3)); 9476528b96dSMatthew G. Knepley PetscFunctionReturn(0); 9486528b96dSMatthew G. Knepley } 9496528b96dSMatthew G. Knepley 95006ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 9516528b96dSMatthew G. Knepley PetscInt n0, 9526528b96dSMatthew G. Knepley void (**g0)(PetscInt, PetscInt, PetscInt, 9536528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9546528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9556528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9566528b96dSMatthew G. Knepley PetscInt n1, 9576528b96dSMatthew G. Knepley void (**g1)(PetscInt, PetscInt, PetscInt, 9586528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9596528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9606528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9616528b96dSMatthew G. Knepley PetscInt n2, 9626528b96dSMatthew G. Knepley void (**g2)(PetscInt, PetscInt, PetscInt, 9636528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9646528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9656528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9666528b96dSMatthew G. Knepley PetscInt n3, 9676528b96dSMatthew G. Knepley void (**g3)(PetscInt, PetscInt, PetscInt, 9686528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9696528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9706528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 9716528b96dSMatthew G. Knepley { 9726528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 9736528b96dSMatthew G. Knepley 9746528b96dSMatthew G. Knepley PetscFunctionBegin; 9759566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (**)(void)) g0)); 9769566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (**)(void)) g1)); 9779566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (**)(void)) g2)); 9789566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (**)(void)) g3)); 9796528b96dSMatthew G. Knepley PetscFunctionReturn(0); 9806528b96dSMatthew G. Knepley } 9816528b96dSMatthew G. Knepley 98206ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetIndexBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 9836528b96dSMatthew G. Knepley PetscInt i0, 9846528b96dSMatthew G. Knepley void (*g0)(PetscInt, PetscInt, PetscInt, 9856528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9866528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9876528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9886528b96dSMatthew G. Knepley PetscInt i1, 9896528b96dSMatthew G. Knepley void (*g1)(PetscInt, PetscInt, PetscInt, 9906528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9916528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9926528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9936528b96dSMatthew G. Knepley PetscInt i2, 9946528b96dSMatthew G. Knepley void (*g2)(PetscInt, PetscInt, PetscInt, 9956528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9966528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9976528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9986528b96dSMatthew G. Knepley PetscInt i3, 9996528b96dSMatthew G. Knepley void (*g3)(PetscInt, PetscInt, PetscInt, 10006528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10016528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10026528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 10036528b96dSMatthew G. Knepley { 10046528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 10056528b96dSMatthew G. Knepley 10066528b96dSMatthew G. Knepley PetscFunctionBegin; 10079566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, i0, (void (*)(void)) g0)); 10089566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, i1, (void (*)(void)) g1)); 10099566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, i2, (void (*)(void)) g2)); 10109566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, i3, (void (*)(void)) g3)); 10116528b96dSMatthew G. Knepley PetscFunctionReturn(0); 10126528b96dSMatthew G. Knepley } 10136528b96dSMatthew G. Knepley 10146528b96dSMatthew G. Knepley PetscErrorCode PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre) 10156528b96dSMatthew G. Knepley { 10166528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 10176528b96dSMatthew G. Knepley 10186528b96dSMatthew G. Knepley PetscFunctionBegin; 10196528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 10206528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJacPre, 2); 10219566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP0], &n0)); 10229566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP1], &n1)); 10239566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP2], &n2)); 10249566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP3], &n3)); 10256528b96dSMatthew G. Knepley *hasJacPre = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE; 10266528b96dSMatthew G. Knepley PetscFunctionReturn(0); 10276528b96dSMatthew G. Knepley } 10286528b96dSMatthew G. Knepley 102906ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 10306528b96dSMatthew G. Knepley PetscInt *n0, 10316528b96dSMatthew G. Knepley void (***g0)(PetscInt, PetscInt, PetscInt, 10326528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10336528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10346528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 10356528b96dSMatthew G. Knepley PetscInt *n1, 10366528b96dSMatthew G. Knepley void (***g1)(PetscInt, PetscInt, PetscInt, 10376528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10386528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10396528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 10406528b96dSMatthew G. Knepley PetscInt *n2, 10416528b96dSMatthew G. Knepley void (***g2)(PetscInt, PetscInt, PetscInt, 10426528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10436528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10446528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 10456528b96dSMatthew G. Knepley PetscInt *n3, 10466528b96dSMatthew G. Knepley void (***g3)(PetscInt, PetscInt, PetscInt, 10476528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10486528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10496528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 10506528b96dSMatthew G. Knepley { 10516528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 10526528b96dSMatthew G. Knepley 10536528b96dSMatthew G. Knepley PetscFunctionBegin; 10549566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (***)(void)) g0)); 10559566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (***)(void)) g1)); 10569566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (***)(void)) g2)); 10579566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (***)(void)) g3)); 10586528b96dSMatthew G. Knepley PetscFunctionReturn(0); 10596528b96dSMatthew G. Knepley } 10606528b96dSMatthew G. Knepley 106106ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormAddBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 10626528b96dSMatthew G. Knepley void (*g0)(PetscInt, PetscInt, PetscInt, 10636528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10646528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10656528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 10666528b96dSMatthew G. Knepley void (*g1)(PetscInt, PetscInt, PetscInt, 10676528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10686528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10696528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 10706528b96dSMatthew G. Knepley void (*g2)(PetscInt, PetscInt, PetscInt, 10716528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10726528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10736528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 10746528b96dSMatthew G. Knepley void (*g3)(PetscInt, PetscInt, PetscInt, 10756528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10766528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10776528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 10786528b96dSMatthew G. Knepley { 10796528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 10806528b96dSMatthew G. Knepley 10816528b96dSMatthew G. Knepley PetscFunctionBegin; 10829566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, (void (*)(void)) g0)); 10839566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, (void (*)(void)) g1)); 10849566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, (void (*)(void)) g2)); 10859566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, (void (*)(void)) g3)); 10866528b96dSMatthew G. Knepley PetscFunctionReturn(0); 10876528b96dSMatthew G. Knepley } 10886528b96dSMatthew G. Knepley 108906ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 10906528b96dSMatthew G. Knepley PetscInt n0, 10916528b96dSMatthew G. Knepley void (**g0)(PetscInt, PetscInt, PetscInt, 10926528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10936528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10946528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 10956528b96dSMatthew G. Knepley PetscInt n1, 10966528b96dSMatthew G. Knepley void (**g1)(PetscInt, PetscInt, PetscInt, 10976528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10986528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 10996528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 11006528b96dSMatthew G. Knepley PetscInt n2, 11016528b96dSMatthew G. Knepley void (**g2)(PetscInt, PetscInt, PetscInt, 11026528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11036528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11046528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 11056528b96dSMatthew G. Knepley PetscInt n3, 11066528b96dSMatthew G. Knepley void (**g3)(PetscInt, PetscInt, PetscInt, 11076528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11086528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11096528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 11106528b96dSMatthew G. Knepley { 11116528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 11126528b96dSMatthew G. Knepley 11136528b96dSMatthew G. Knepley PetscFunctionBegin; 11149566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (**)(void)) g0)); 11159566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (**)(void)) g1)); 11169566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (**)(void)) g2)); 11179566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (**)(void)) g3)); 11186528b96dSMatthew G. Knepley PetscFunctionReturn(0); 11196528b96dSMatthew G. Knepley } 11206528b96dSMatthew G. Knepley 112106ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetIndexBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 11226528b96dSMatthew G. Knepley PetscInt i0, 11236528b96dSMatthew G. Knepley void (*g0)(PetscInt, PetscInt, PetscInt, 11246528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11256528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11266528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 11276528b96dSMatthew G. Knepley PetscInt i1, 11286528b96dSMatthew G. Knepley void (*g1)(PetscInt, PetscInt, PetscInt, 11296528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11306528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11316528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 11326528b96dSMatthew G. Knepley PetscInt i2, 11336528b96dSMatthew G. Knepley void (*g2)(PetscInt, PetscInt, PetscInt, 11346528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11356528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11366528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 11376528b96dSMatthew G. Knepley PetscInt i3, 11386528b96dSMatthew G. Knepley void (*g3)(PetscInt, PetscInt, PetscInt, 11396528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11406528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11416528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 11426528b96dSMatthew G. Knepley { 11436528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 11446528b96dSMatthew G. Knepley 11456528b96dSMatthew G. Knepley PetscFunctionBegin; 11469566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, i0, (void (*)(void)) g0)); 11479566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, i1, (void (*)(void)) g1)); 11489566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, i2, (void (*)(void)) g2)); 11499566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, i3, (void (*)(void)) g3)); 11506528b96dSMatthew G. Knepley PetscFunctionReturn(0); 11516528b96dSMatthew G. Knepley } 11526528b96dSMatthew G. Knepley 11536528b96dSMatthew G. Knepley PetscErrorCode PetscWeakFormHasDynamicJacobian(PetscWeakForm wf, PetscBool *hasDynJac) 11546528b96dSMatthew G. Knepley { 11556528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 11566528b96dSMatthew G. Knepley 11576528b96dSMatthew G. Knepley PetscFunctionBegin; 11586528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 11596528b96dSMatthew G. Knepley PetscValidBoolPointer(hasDynJac, 2); 11609566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT0], &n0)); 11619566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT1], &n1)); 11629566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT2], &n2)); 11639566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT3], &n3)); 11646528b96dSMatthew G. Knepley *hasDynJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE; 11656528b96dSMatthew G. Knepley PetscFunctionReturn(0); 11666528b96dSMatthew G. Knepley } 11676528b96dSMatthew G. Knepley 116806ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 11696528b96dSMatthew G. Knepley PetscInt *n0, 11706528b96dSMatthew G. Knepley void (***g0)(PetscInt, PetscInt, PetscInt, 11716528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11726528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11736528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 11746528b96dSMatthew G. Knepley PetscInt *n1, 11756528b96dSMatthew G. Knepley void (***g1)(PetscInt, PetscInt, PetscInt, 11766528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11776528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11786528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 11796528b96dSMatthew G. Knepley PetscInt *n2, 11806528b96dSMatthew G. Knepley void (***g2)(PetscInt, PetscInt, PetscInt, 11816528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11826528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11836528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 11846528b96dSMatthew G. Knepley PetscInt *n3, 11856528b96dSMatthew G. Knepley void (***g3)(PetscInt, PetscInt, PetscInt, 11866528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11876528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 11886528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 11896528b96dSMatthew G. Knepley { 11906528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 11916528b96dSMatthew G. Knepley 11926528b96dSMatthew G. Knepley PetscFunctionBegin; 11939566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (***)(void)) g0)); 11949566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (***)(void)) g1)); 11959566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (***)(void)) g2)); 11969566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (***)(void)) g3)); 11976528b96dSMatthew G. Knepley PetscFunctionReturn(0); 11986528b96dSMatthew G. Knepley } 11996528b96dSMatthew G. Knepley 120006ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormAddDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 12016528b96dSMatthew G. Knepley void (*g0)(PetscInt, PetscInt, PetscInt, 12026528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12036528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12046528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 12056528b96dSMatthew G. Knepley void (*g1)(PetscInt, PetscInt, PetscInt, 12066528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12076528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12086528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 12096528b96dSMatthew G. Knepley void (*g2)(PetscInt, PetscInt, PetscInt, 12106528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12116528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12126528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 12136528b96dSMatthew G. Knepley void (*g3)(PetscInt, PetscInt, PetscInt, 12146528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12156528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12166528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 12176528b96dSMatthew G. Knepley { 12186528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 12196528b96dSMatthew G. Knepley 12206528b96dSMatthew G. Knepley PetscFunctionBegin; 12219566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, (void (*)(void)) g0)); 12229566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, (void (*)(void)) g1)); 12239566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, (void (*)(void)) g2)); 12249566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, (void (*)(void)) g3)); 12256528b96dSMatthew G. Knepley PetscFunctionReturn(0); 12266528b96dSMatthew G. Knepley } 12276528b96dSMatthew G. Knepley 122806ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 12296528b96dSMatthew G. Knepley PetscInt n0, 12306528b96dSMatthew G. Knepley void (**g0)(PetscInt, PetscInt, PetscInt, 12316528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12326528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12336528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 12346528b96dSMatthew G. Knepley PetscInt n1, 12356528b96dSMatthew G. Knepley void (**g1)(PetscInt, PetscInt, PetscInt, 12366528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12376528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12386528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 12396528b96dSMatthew G. Knepley PetscInt n2, 12406528b96dSMatthew G. Knepley void (**g2)(PetscInt, PetscInt, PetscInt, 12416528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12426528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12436528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 12446528b96dSMatthew G. Knepley PetscInt n3, 12456528b96dSMatthew G. Knepley void (**g3)(PetscInt, PetscInt, PetscInt, 12466528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12476528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12486528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 12496528b96dSMatthew G. Knepley { 12506528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 12516528b96dSMatthew G. Knepley 12526528b96dSMatthew G. Knepley PetscFunctionBegin; 12539566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (**)(void)) g0)); 12549566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (**)(void)) g1)); 12559566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (**)(void)) g2)); 12569566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (**)(void)) g3)); 12576528b96dSMatthew G. Knepley PetscFunctionReturn(0); 12586528b96dSMatthew G. Knepley } 12596528b96dSMatthew G. Knepley 126006ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetIndexDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, 12616528b96dSMatthew G. Knepley PetscInt i0, 12626528b96dSMatthew G. Knepley void (*g0)(PetscInt, PetscInt, PetscInt, 12636528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12646528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12656528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 12666528b96dSMatthew G. Knepley PetscInt i1, 12676528b96dSMatthew G. Knepley void (*g1)(PetscInt, PetscInt, PetscInt, 12686528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12696528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12706528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 12716528b96dSMatthew G. Knepley PetscInt i2, 12726528b96dSMatthew G. Knepley void (*g2)(PetscInt, PetscInt, PetscInt, 12736528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12746528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12756528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 12766528b96dSMatthew G. Knepley PetscInt i3, 12776528b96dSMatthew G. Knepley void (*g3)(PetscInt, PetscInt, PetscInt, 12786528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12796528b96dSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 12806528b96dSMatthew G. Knepley PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 12816528b96dSMatthew G. Knepley { 12826528b96dSMatthew G. Knepley PetscInt find = f*wf->Nf + g; 12836528b96dSMatthew G. Knepley 12846528b96dSMatthew G. Knepley PetscFunctionBegin; 12859566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, i0, (void (*)(void)) g0)); 12869566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, i1, (void (*)(void)) g1)); 12879566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, i2, (void (*)(void)) g2)); 12889566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, i3, (void (*)(void)) g3)); 12896528b96dSMatthew G. Knepley PetscFunctionReturn(0); 12906528b96dSMatthew G. Knepley } 12916528b96dSMatthew G. Knepley 129206ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormGetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n, 12936528b96dSMatthew G. Knepley void (***r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)) 12946528b96dSMatthew G. Knepley { 12956528b96dSMatthew G. Knepley PetscFunctionBegin; 12969566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (***)(void)) r)); 12976528b96dSMatthew G. Knepley PetscFunctionReturn(0); 12986528b96dSMatthew G. Knepley } 12996528b96dSMatthew G. Knepley 130006ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 13016528b96dSMatthew G. Knepley PetscInt n, 13026528b96dSMatthew G. Knepley void (**r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)) 13036528b96dSMatthew G. Knepley { 13046528b96dSMatthew G. Knepley PetscFunctionBegin; 13059566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (**)(void)) r)); 13066528b96dSMatthew G. Knepley PetscFunctionReturn(0); 13076528b96dSMatthew G. Knepley } 13086528b96dSMatthew G. Knepley 130906ad1575SMatthew G. Knepley PetscErrorCode PetscWeakFormSetIndexRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, 13106528b96dSMatthew G. Knepley PetscInt i, 13116528b96dSMatthew G. Knepley void (*r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)) 13126528b96dSMatthew G. Knepley { 13136528b96dSMatthew G. Knepley PetscFunctionBegin; 13149566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, i, (void (*)(void)) r)); 13156528b96dSMatthew G. Knepley PetscFunctionReturn(0); 13166528b96dSMatthew G. Knepley } 13176528b96dSMatthew G. Knepley 13186528b96dSMatthew G. Knepley /*@ 13196528b96dSMatthew G. Knepley PetscWeakFormGetNumFields - Returns the number of fields 13206528b96dSMatthew G. Knepley 13216528b96dSMatthew G. Knepley Not collective 13226528b96dSMatthew G. Knepley 13236528b96dSMatthew G. Knepley Input Parameter: 13246528b96dSMatthew G. Knepley . wf - The PetscWeakForm object 13256528b96dSMatthew G. Knepley 13266528b96dSMatthew G. Knepley Output Parameter: 1327a5b23f4aSJose E. Roman . Nf - The number of fields 13286528b96dSMatthew G. Knepley 13296528b96dSMatthew G. Knepley Level: beginner 13306528b96dSMatthew G. Knepley 1331db781477SPatrick Sanan .seealso: `PetscWeakFormSetNumFields()`, `PetscWeakFormCreate()` 13326528b96dSMatthew G. Knepley @*/ 13336528b96dSMatthew G. Knepley PetscErrorCode PetscWeakFormGetNumFields(PetscWeakForm wf, PetscInt *Nf) 13346528b96dSMatthew G. Knepley { 13356528b96dSMatthew G. Knepley PetscFunctionBegin; 13366528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 1337dadcf809SJacob Faibussowitsch PetscValidIntPointer(Nf, 2); 13386528b96dSMatthew G. Knepley *Nf = wf->Nf; 13396528b96dSMatthew G. Knepley PetscFunctionReturn(0); 13406528b96dSMatthew G. Knepley } 13416528b96dSMatthew G. Knepley 13426528b96dSMatthew G. Knepley /*@ 13436528b96dSMatthew G. Knepley PetscWeakFormSetNumFields - Sets the number of fields 13446528b96dSMatthew G. Knepley 13456528b96dSMatthew G. Knepley Not collective 13466528b96dSMatthew G. Knepley 13476528b96dSMatthew G. Knepley Input Parameters: 13486528b96dSMatthew G. Knepley + wf - The PetscWeakForm object 13496528b96dSMatthew G. Knepley - Nf - The number of fields 13506528b96dSMatthew G. Knepley 13516528b96dSMatthew G. Knepley Level: beginner 13526528b96dSMatthew G. Knepley 1353db781477SPatrick Sanan .seealso: `PetscWeakFormGetNumFields()`, `PetscWeakFormCreate()` 13546528b96dSMatthew G. Knepley @*/ 13556528b96dSMatthew G. Knepley PetscErrorCode PetscWeakFormSetNumFields(PetscWeakForm wf, PetscInt Nf) 13566528b96dSMatthew G. Knepley { 13576528b96dSMatthew G. Knepley PetscFunctionBegin; 13586528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 13596528b96dSMatthew G. Knepley wf->Nf = Nf; 13606528b96dSMatthew G. Knepley PetscFunctionReturn(0); 13616528b96dSMatthew G. Knepley } 13626528b96dSMatthew G. Knepley 13636528b96dSMatthew G. Knepley /*@ 13646528b96dSMatthew G. Knepley PetscWeakFormDestroy - Destroys a PetscWeakForm object 13656528b96dSMatthew G. Knepley 13666528b96dSMatthew G. Knepley Collective on wf 13676528b96dSMatthew G. Knepley 13686528b96dSMatthew G. Knepley Input Parameter: 13696528b96dSMatthew G. Knepley . wf - the PetscWeakForm object to destroy 13706528b96dSMatthew G. Knepley 13716528b96dSMatthew G. Knepley Level: developer 13726528b96dSMatthew G. Knepley 1373db781477SPatrick Sanan .seealso `PetscWeakFormCreate()`, `PetscWeakFormView()` 13746528b96dSMatthew G. Knepley @*/ 13756528b96dSMatthew G. Knepley PetscErrorCode PetscWeakFormDestroy(PetscWeakForm *wf) 13766528b96dSMatthew G. Knepley { 137706ad1575SMatthew G. Knepley PetscInt f; 13786528b96dSMatthew G. Knepley 13796528b96dSMatthew G. Knepley PetscFunctionBegin; 13806528b96dSMatthew G. Knepley if (!*wf) PetscFunctionReturn(0); 13816528b96dSMatthew G. Knepley PetscValidHeaderSpecific((*wf), PETSCWEAKFORM_CLASSID, 1); 13826528b96dSMatthew G. Knepley 13836528b96dSMatthew G. Knepley if (--((PetscObject)(*wf))->refct > 0) {*wf = NULL; PetscFunctionReturn(0);} 13846528b96dSMatthew G. Knepley ((PetscObject) (*wf))->refct = 0; 13859566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferDestroy(&(*wf)->funcs)); 13869566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormDestroy(&(*wf)->form[f])); 13879566063dSJacob Faibussowitsch PetscCall(PetscFree((*wf)->form)); 13889566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(wf)); 13896528b96dSMatthew G. Knepley PetscFunctionReturn(0); 13906528b96dSMatthew G. Knepley } 13916528b96dSMatthew G. Knepley 139245480ffeSMatthew G. Knepley static PetscErrorCode PetscWeakFormViewTable_Ascii(PetscWeakForm wf, PetscViewer viewer, PetscBool splitField, const char tableName[], PetscHMapForm map) 13936528b96dSMatthew G. Knepley { 139445480ffeSMatthew G. Knepley PetscInt Nf = wf->Nf, Nk, k; 13956528b96dSMatthew G. Knepley 13966528b96dSMatthew G. Knepley PetscFunctionBegin; 13979566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(map, &Nk)); 13986528b96dSMatthew G. Knepley if (Nk) { 139906ad1575SMatthew G. Knepley PetscFormKey *keys; 1400165f9cc3SJed Brown void (**funcs)(void) = NULL; 14015fedec97SMatthew G. Knepley const char **names; 14025fedec97SMatthew G. Knepley PetscInt *values, *idx1, *idx2, *idx; 14035fedec97SMatthew G. Knepley PetscBool showPart = PETSC_FALSE, showPointer = PETSC_FALSE; 14045fedec97SMatthew G. Knepley PetscInt off = 0; 14056528b96dSMatthew G. Knepley 14069566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(Nk, &keys, Nk, &names, Nk, &values, Nk, &idx1, Nk, &idx2, Nk, &idx)); 14079566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetKeys(map, &off, keys)); 14085fedec97SMatthew G. Knepley /* Sort keys by label name and value */ 14095fedec97SMatthew G. Knepley { 14105fedec97SMatthew G. Knepley /* First sort values */ 14115fedec97SMatthew G. Knepley for (k = 0; k < Nk; ++k) {values[k] = keys[k].value; idx1[k] = k;} 14129566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(Nk, values, idx1)); 14135fedec97SMatthew G. Knepley /* If the string sort is stable, it will be sorted correctly overall */ 14145fedec97SMatthew G. Knepley for (k = 0; k < Nk; ++k) { 14159566063dSJacob Faibussowitsch if (keys[idx1[k]].label) PetscCall(PetscObjectGetName((PetscObject) keys[idx1[k]].label, &names[k])); 14165fedec97SMatthew G. Knepley else {names[k] = "";} 14175fedec97SMatthew G. Knepley idx2[k] = k; 14185fedec97SMatthew G. Knepley } 14199566063dSJacob Faibussowitsch PetscCall(PetscSortStrWithPermutation(Nk, names, idx2)); 14205fedec97SMatthew G. Knepley for (k = 0; k < Nk; ++k) { 14219566063dSJacob Faibussowitsch if (keys[k].label) PetscCall(PetscObjectGetName((PetscObject) keys[k].label, &names[k])); 14225fedec97SMatthew G. Knepley else {names[k] = "";} 14235fedec97SMatthew G. Knepley idx[k] = idx1[idx2[k]]; 14245fedec97SMatthew G. Knepley } 14255fedec97SMatthew G. Knepley } 14269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", tableName)); 14279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 14286528b96dSMatthew G. Knepley for (k = 0; k < Nk; ++k) { 14290944a4d6SMatthew G. Knepley if (keys[k].part != 0) showPart = PETSC_TRUE; 14300944a4d6SMatthew G. Knepley } 14310944a4d6SMatthew G. Knepley for (k = 0; k < Nk; ++k) { 14325fedec97SMatthew G. Knepley const PetscInt i = idx[k]; 14335fedec97SMatthew G. Knepley PetscInt n, f; 14345fedec97SMatthew G. Knepley 14355fedec97SMatthew G. Knepley if (keys[i].label) { 143663a3b9bcSJacob Faibussowitsch if (showPointer) PetscCall(PetscViewerASCIIPrintf(viewer, "(%s:%p, %" PetscInt_FMT ") ", names[i], keys[i].label, keys[i].value)); 143763a3b9bcSJacob Faibussowitsch else PetscCall(PetscViewerASCIIPrintf(viewer, "(%s, %" PetscInt_FMT ") ", names[i], keys[i].value)); 143863a3b9bcSJacob Faibussowitsch } 14399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 144063a3b9bcSJacob Faibussowitsch if (splitField) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ", %" PetscInt_FMT ") ", keys[i].field/Nf, keys[i].field%Nf)); 144163a3b9bcSJacob Faibussowitsch else PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") ", keys[i].field)); 144263a3b9bcSJacob Faibussowitsch if (showPart) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") ", keys[i].part)); 14439566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, map, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &n, &funcs)); 14445fedec97SMatthew G. Knepley for (f = 0; f < n; ++f) { 1445258ec3d2SMatthew G. Knepley char *fname; 14465fedec97SMatthew G. Knepley size_t len, l; 1447258ec3d2SMatthew G. Knepley 14489566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 14499566063dSJacob Faibussowitsch PetscCall(PetscDLAddr(funcs[f], &fname)); 14505fedec97SMatthew G. Knepley if (fname) { 14515fedec97SMatthew G. Knepley /* Eliminate argument types */ 14529566063dSJacob Faibussowitsch PetscCall(PetscStrlen(fname, &len)); 14535fedec97SMatthew G. Knepley for (l = 0; l < len; ++l) if (fname[l] == '(') {fname[l] = '\0'; break;} 14549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%s", fname)); 14555fedec97SMatthew G. Knepley } else if (showPointer) { 14569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%p", funcs[f])); 14575fedec97SMatthew G. Knepley } 14589566063dSJacob Faibussowitsch PetscCall(PetscFree(fname)); 14596528b96dSMatthew G. Knepley } 14609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 14619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 14626528b96dSMatthew G. Knepley } 14639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 14649566063dSJacob Faibussowitsch PetscCall(PetscFree6(keys, names, values, idx1, idx2, idx)); 14656528b96dSMatthew G. Knepley } 14666528b96dSMatthew G. Knepley PetscFunctionReturn(0); 14676528b96dSMatthew G. Knepley } 14686528b96dSMatthew G. Knepley 14696528b96dSMatthew G. Knepley static PetscErrorCode PetscWeakFormView_Ascii(PetscWeakForm wf, PetscViewer viewer) 14706528b96dSMatthew G. Knepley { 14716528b96dSMatthew G. Knepley PetscViewerFormat format; 147206ad1575SMatthew G. Knepley PetscInt f; 14736528b96dSMatthew G. Knepley 14746528b96dSMatthew G. Knepley PetscFunctionBegin; 14759566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 147663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Weak Form System with %" PetscInt_FMT " fields\n", wf->Nf)); 14779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 147806ad1575SMatthew G. Knepley for (f = 0; f < PETSC_NUM_WF; ++f) { 14799566063dSJacob Faibussowitsch PetscCall(PetscWeakFormViewTable_Ascii(wf, viewer, PETSC_TRUE, PetscWeakFormKinds[f], wf->form[f])); 148006ad1575SMatthew G. Knepley } 14819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 14826528b96dSMatthew G. Knepley PetscFunctionReturn(0); 14836528b96dSMatthew G. Knepley } 14846528b96dSMatthew G. Knepley 14856528b96dSMatthew G. Knepley /*@C 14866528b96dSMatthew G. Knepley PetscWeakFormView - Views a PetscWeakForm 14876528b96dSMatthew G. Knepley 14886528b96dSMatthew G. Knepley Collective on wf 14896528b96dSMatthew G. Knepley 1490d8d19677SJose E. Roman Input Parameters: 14916528b96dSMatthew G. Knepley + wf - the PetscWeakForm object to view 14926528b96dSMatthew G. Knepley - v - the viewer 14936528b96dSMatthew G. Knepley 14946528b96dSMatthew G. Knepley Level: developer 14956528b96dSMatthew G. Knepley 1496db781477SPatrick Sanan .seealso `PetscWeakFormDestroy()`, `PetscWeakFormCreate()` 14976528b96dSMatthew G. Knepley @*/ 14986528b96dSMatthew G. Knepley PetscErrorCode PetscWeakFormView(PetscWeakForm wf, PetscViewer v) 14996528b96dSMatthew G. Knepley { 15006528b96dSMatthew G. Knepley PetscBool iascii; 15016528b96dSMatthew G. Knepley 15026528b96dSMatthew G. Knepley PetscFunctionBegin; 15036528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 15049566063dSJacob Faibussowitsch if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) wf), &v)); 15056528b96dSMatthew G. Knepley else {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);} 15069566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii)); 15079566063dSJacob Faibussowitsch if (iascii) PetscCall(PetscWeakFormView_Ascii(wf, v)); 1508*dbbe0bcdSBarry Smith PetscTryTypeMethod(wf,view, v); 15096528b96dSMatthew G. Knepley PetscFunctionReturn(0); 15106528b96dSMatthew G. Knepley } 15116528b96dSMatthew G. Knepley 15126528b96dSMatthew G. Knepley /*@ 15136528b96dSMatthew G. Knepley PetscWeakFormCreate - Creates an empty PetscWeakForm object. 15146528b96dSMatthew G. Knepley 15156528b96dSMatthew G. Knepley Collective 15166528b96dSMatthew G. Knepley 15176528b96dSMatthew G. Knepley Input Parameter: 15186528b96dSMatthew G. Knepley . comm - The communicator for the PetscWeakForm object 15196528b96dSMatthew G. Knepley 15206528b96dSMatthew G. Knepley Output Parameter: 15216528b96dSMatthew G. Knepley . wf - The PetscWeakForm object 15226528b96dSMatthew G. Knepley 15236528b96dSMatthew G. Knepley Level: beginner 15246528b96dSMatthew G. Knepley 1525db781477SPatrick Sanan .seealso: `PetscDS`, `PetscWeakFormDestroy()` 15266528b96dSMatthew G. Knepley @*/ 15276528b96dSMatthew G. Knepley PetscErrorCode PetscWeakFormCreate(MPI_Comm comm, PetscWeakForm *wf) 15286528b96dSMatthew G. Knepley { 15296528b96dSMatthew G. Knepley PetscWeakForm p; 153006ad1575SMatthew G. Knepley PetscInt f; 15316528b96dSMatthew G. Knepley 15326528b96dSMatthew G. Knepley PetscFunctionBegin; 15336528b96dSMatthew G. Knepley PetscValidPointer(wf, 2); 15346528b96dSMatthew G. Knepley *wf = NULL; 15359566063dSJacob Faibussowitsch PetscCall(PetscDSInitializePackage()); 15366528b96dSMatthew G. Knepley 15379566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(p, PETSCWEAKFORM_CLASSID, "PetscWeakForm", "Weak Form System", "PetscWeakForm", comm, PetscWeakFormDestroy, PetscWeakFormView)); 15386528b96dSMatthew G. Knepley 15396528b96dSMatthew G. Knepley p->Nf = 0; 15409566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreate(sizeof(&PetscWeakFormCreate), 2, &p->funcs)); 15419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PETSC_NUM_WF, &p->form)); 15429566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormCreate(&p->form[f])); 15436528b96dSMatthew G. Knepley *wf = p; 15446528b96dSMatthew G. Knepley PetscFunctionReturn(0); 15456528b96dSMatthew G. Knepley } 1546