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 7d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscChunkBufferCreate(size_t unitbytes, size_t expected, PetscChunkBuffer **buffer) 8d71ae5a4SJacob Faibussowitsch { 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 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscChunkBufferDuplicate(PetscChunkBuffer *buffer, PetscChunkBuffer **bufferNew) 19d71ae5a4SJacob Faibussowitsch { 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 30d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscChunkBufferDestroy(PetscChunkBuffer **buffer) 31d71ae5a4SJacob Faibussowitsch { 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 38d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscChunkBufferCreateChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk) 39d71ae5a4SJacob Faibussowitsch { 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 58d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscChunkBufferEnlargeChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk) 59d71ae5a4SJacob Faibussowitsch { 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 @*/ 92d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFormKeySort(PetscInt n, PetscFormKey arr[]) 93d71ae5a4SJacob Faibussowitsch { 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 101d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt *n, void (***func)()) 102d71ae5a4SJacob Faibussowitsch { 10306ad1575SMatthew G. Knepley PetscFormKey key; 1046528b96dSMatthew G. Knepley PetscChunk chunk; 1056528b96dSMatthew G. Knepley 1066528b96dSMatthew G. Knepley PetscFunctionBegin; 1079371c9d4SSatish Balay key.label = label; 1089371c9d4SSatish Balay key.value = value; 1099371c9d4SSatish Balay key.field = f; 1109371c9d4SSatish Balay key.part = part; 1119566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1129371c9d4SSatish Balay if (chunk.size < 0) { 1139371c9d4SSatish Balay *n = 0; 1149371c9d4SSatish Balay *func = NULL; 1159371c9d4SSatish Balay } else { 1169371c9d4SSatish Balay *n = chunk.size; 1179371c9d4SSatish Balay *func = (void (**)()) & wf->funcs->array[chunk.start]; 1189371c9d4SSatish Balay } 1196528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1206528b96dSMatthew G. Knepley } 1216528b96dSMatthew G. Knepley 1226528b96dSMatthew G. Knepley /* A NULL argument for func causes this to clear the key */ 123d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt n, void (**func)()) 124d71ae5a4SJacob Faibussowitsch { 12506ad1575SMatthew G. Knepley PetscFormKey key; 1266528b96dSMatthew G. Knepley PetscChunk chunk; 1276528b96dSMatthew G. Knepley PetscInt i; 1286528b96dSMatthew G. Knepley 1296528b96dSMatthew G. Knepley PetscFunctionBegin; 1309371c9d4SSatish Balay key.label = label; 1319371c9d4SSatish Balay key.value = value; 1329371c9d4SSatish Balay key.field = f; 1339371c9d4SSatish Balay key.part = part; 1346528b96dSMatthew G. Knepley if (!func) { 1359566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDel(ht, key)); 1366528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1371baa6e33SBarry Smith } else PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1386528b96dSMatthew G. Knepley if (chunk.size < 0) { 1399566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(wf->funcs, n, &chunk)); 1409566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1416528b96dSMatthew G. Knepley } else if (chunk.size <= n) { 1429566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, n - chunk.size, &chunk)); 1439566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1446528b96dSMatthew G. Knepley } 1452baeeceeSMatthew G. Knepley for (i = 0; i < n; ++i) ((void (**)()) & wf->funcs->array[chunk.start])[i] = func[i]; 1466528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1476528b96dSMatthew G. Knepley } 1486528b96dSMatthew G. Knepley 149d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormAddFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, void (*func)()) 150d71ae5a4SJacob Faibussowitsch { 15106ad1575SMatthew G. Knepley PetscFormKey key; 1526528b96dSMatthew G. Knepley PetscChunk chunk; 1536528b96dSMatthew G. Knepley 1546528b96dSMatthew G. Knepley PetscFunctionBegin; 1556528b96dSMatthew G. Knepley if (!func) PetscFunctionReturn(0); 1569371c9d4SSatish Balay key.label = label; 1579371c9d4SSatish Balay key.value = value; 1589371c9d4SSatish Balay key.field = f; 1599371c9d4SSatish Balay key.part = part; 1609566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1616528b96dSMatthew G. Knepley if (chunk.size < 0) { 1629566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(wf->funcs, 1, &chunk)); 1639566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1642baeeceeSMatthew G. Knepley ((void (**)()) & wf->funcs->array[chunk.start])[0] = func; 1656528b96dSMatthew G. Knepley } else { 1669566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, 1, &chunk)); 1679566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1682baeeceeSMatthew G. Knepley ((void (**)()) & wf->funcs->array[chunk.start])[chunk.size - 1] = func; 1696528b96dSMatthew G. Knepley } 1706528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1716528b96dSMatthew G. Knepley } 1726528b96dSMatthew G. Knepley 173d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, void (**func)()) 174d71ae5a4SJacob Faibussowitsch { 17506ad1575SMatthew G. Knepley PetscFormKey key; 1766528b96dSMatthew G. Knepley PetscChunk chunk; 1776528b96dSMatthew G. Knepley 1786528b96dSMatthew G. Knepley PetscFunctionBegin; 1799371c9d4SSatish Balay key.label = label; 1809371c9d4SSatish Balay key.value = value; 1819371c9d4SSatish Balay key.field = f; 1829371c9d4SSatish Balay key.part = part; 1839566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1849371c9d4SSatish Balay if (chunk.size < 0) { 1859371c9d4SSatish Balay *func = NULL; 1869371c9d4SSatish Balay } else { 18763a3b9bcSJacob Faibussowitsch PetscCheck(ind < chunk.size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", ind, chunk.size); 1882baeeceeSMatthew G. Knepley *func = ((void (**)()) & wf->funcs->array[chunk.start])[ind]; 1896528b96dSMatthew G. Knepley } 1906528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1916528b96dSMatthew G. Knepley } 1926528b96dSMatthew G. Knepley 19306ad1575SMatthew G. Knepley /* Ignore a NULL func */ 194d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, void (*func)()) 195d71ae5a4SJacob Faibussowitsch { 19606ad1575SMatthew G. Knepley PetscFormKey key; 1976528b96dSMatthew G. Knepley PetscChunk chunk; 1986528b96dSMatthew G. Knepley 1996528b96dSMatthew G. Knepley PetscFunctionBegin; 20006ad1575SMatthew G. Knepley if (!func) PetscFunctionReturn(0); 2019371c9d4SSatish Balay key.label = label; 2029371c9d4SSatish Balay key.value = value; 2039371c9d4SSatish Balay key.field = f; 2049371c9d4SSatish Balay key.part = part; 2059566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 2066528b96dSMatthew G. Knepley if (chunk.size < 0) { 2079566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(wf->funcs, ind + 1, &chunk)); 2089566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 2096528b96dSMatthew G. Knepley } else if (chunk.size <= ind) { 2109566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, ind - chunk.size + 1, &chunk)); 2119566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 2126528b96dSMatthew G. Knepley } 2132baeeceeSMatthew G. Knepley ((void (**)()) & wf->funcs->array[chunk.start])[ind] = func; 2146528b96dSMatthew G. Knepley PetscFunctionReturn(0); 2156528b96dSMatthew G. Knepley } 2166528b96dSMatthew G. Knepley 217d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormClearIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind) 218d71ae5a4SJacob Faibussowitsch { 21906ad1575SMatthew G. Knepley PetscFormKey key; 22006ad1575SMatthew G. Knepley PetscChunk chunk; 22106ad1575SMatthew G. Knepley 22206ad1575SMatthew G. Knepley PetscFunctionBegin; 2239371c9d4SSatish Balay key.label = label; 2249371c9d4SSatish Balay key.value = value; 2259371c9d4SSatish Balay key.field = f; 2269371c9d4SSatish Balay key.part = part; 2279566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 22806ad1575SMatthew G. Knepley if (chunk.size < 0) { 22906ad1575SMatthew G. Knepley PetscFunctionReturn(0); 23006ad1575SMatthew G. Knepley } else if (!ind && chunk.size == 1) { 2319566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDel(ht, key)); 23206ad1575SMatthew G. Knepley PetscFunctionReturn(0); 23306ad1575SMatthew G. Knepley } else if (chunk.size <= ind) { 23406ad1575SMatthew G. Knepley PetscFunctionReturn(0); 23506ad1575SMatthew G. Knepley } 23606ad1575SMatthew G. Knepley ((void (**)()) & wf->funcs->array[chunk.start])[ind] = NULL; 23706ad1575SMatthew G. Knepley PetscFunctionReturn(0); 23806ad1575SMatthew G. Knepley } 23906ad1575SMatthew G. Knepley 24045480ffeSMatthew G. Knepley /*@ 24145480ffeSMatthew G. Knepley PetscWeakFormCopy - Copy the pointwise functions to another PetscWeakForm 24245480ffeSMatthew G. Knepley 24345480ffeSMatthew G. Knepley Not Collective 24445480ffeSMatthew G. Knepley 24545480ffeSMatthew G. Knepley Input Parameter: 24645480ffeSMatthew G. Knepley . wf - The original PetscWeakForm 24745480ffeSMatthew G. Knepley 24845480ffeSMatthew G. Knepley Output Parameter: 24945480ffeSMatthew G. Knepley . wfNew - The copy PetscWeakForm 25045480ffeSMatthew G. Knepley 25145480ffeSMatthew G. Knepley Level: intermediate 25245480ffeSMatthew G. Knepley 253db781477SPatrick Sanan .seealso: `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 25445480ffeSMatthew G. Knepley @*/ 255d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormCopy(PetscWeakForm wf, PetscWeakForm wfNew) 256d71ae5a4SJacob Faibussowitsch { 25706ad1575SMatthew G. Knepley PetscInt f; 25845480ffeSMatthew G. Knepley 25945480ffeSMatthew G. Knepley PetscFunctionBegin; 26045480ffeSMatthew G. Knepley wfNew->Nf = wf->Nf; 2619566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferDestroy(&wfNew->funcs)); 2629566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferDuplicate(wf->funcs, &wfNew->funcs)); 26306ad1575SMatthew G. Knepley for (f = 0; f < PETSC_NUM_WF; ++f) { 2649566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDestroy(&wfNew->form[f])); 2659566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDuplicate(wf->form[f], &wfNew->form[f])); 26606ad1575SMatthew G. Knepley } 26706ad1575SMatthew G. Knepley PetscFunctionReturn(0); 26806ad1575SMatthew G. Knepley } 26906ad1575SMatthew G. Knepley 27006ad1575SMatthew G. Knepley /*@ 27106ad1575SMatthew G. Knepley PetscWeakFormClear - Clear all functions from the PetscWeakForm 27206ad1575SMatthew G. Knepley 27306ad1575SMatthew G. Knepley Not Collective 27406ad1575SMatthew G. Knepley 27506ad1575SMatthew G. Knepley Input Parameter: 27606ad1575SMatthew G. Knepley . wf - The original PetscWeakForm 27706ad1575SMatthew G. Knepley 27806ad1575SMatthew G. Knepley Level: intermediate 27906ad1575SMatthew G. Knepley 280db781477SPatrick Sanan .seealso: `PetscWeakFormCopy()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 28106ad1575SMatthew G. Knepley @*/ 282d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormClear(PetscWeakForm wf) 283d71ae5a4SJacob Faibussowitsch { 28406ad1575SMatthew G. Knepley PetscInt f; 28506ad1575SMatthew G. Knepley 28606ad1575SMatthew G. Knepley PetscFunctionBegin; 2879566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormClear(wf->form[f])); 28845480ffeSMatthew G. Knepley PetscFunctionReturn(0); 28945480ffeSMatthew G. Knepley } 29045480ffeSMatthew G. Knepley 291d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscWeakFormRewriteKeys_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label, PetscInt Nv, const PetscInt values[]) 292d71ae5a4SJacob Faibussowitsch { 29306ad1575SMatthew G. Knepley PetscFormKey *keys; 29445480ffeSMatthew G. Knepley PetscInt n, i, v, off = 0; 29545480ffeSMatthew G. Knepley 29645480ffeSMatthew G. Knepley PetscFunctionBegin; 2979566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(hmap, &n)); 2989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &keys)); 2999566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetKeys(hmap, &off, keys)); 30045480ffeSMatthew G. Knepley for (i = 0; i < n; ++i) { 30145480ffeSMatthew G. Knepley if (keys[i].label == label) { 30245480ffeSMatthew G. Knepley PetscBool clear = PETSC_TRUE; 30345480ffeSMatthew G. Knepley void (**funcs)(); 30445480ffeSMatthew G. Knepley PetscInt Nf; 30545480ffeSMatthew G. Knepley 3069566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 30745480ffeSMatthew G. Knepley for (v = 0; v < Nv; ++v) { 3089566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, values[v], keys[i].field, keys[i].part, Nf, funcs)); 30945480ffeSMatthew G. Knepley if (values[v] == keys[i].value) clear = PETSC_FALSE; 31045480ffeSMatthew G. Knepley } 3119566063dSJacob Faibussowitsch if (clear) PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL)); 31245480ffeSMatthew G. Knepley } 31345480ffeSMatthew G. Knepley } 3149566063dSJacob Faibussowitsch PetscCall(PetscFree(keys)); 31545480ffeSMatthew G. Knepley PetscFunctionReturn(0); 31645480ffeSMatthew G. Knepley } 31745480ffeSMatthew G. Knepley 31845480ffeSMatthew G. Knepley /*@C 31945480ffeSMatthew G. Knepley PetscWeakFormRewriteKeys - Change any key on the given label to use the new set of label values 32045480ffeSMatthew G. Knepley 32145480ffeSMatthew G. Knepley Not Collective 32245480ffeSMatthew G. Knepley 32345480ffeSMatthew G. Knepley Input Parameters: 32445480ffeSMatthew G. Knepley + wf - The original PetscWeakForm 32545480ffeSMatthew G. Knepley . label - The label to change keys for 32645480ffeSMatthew G. Knepley . Nv - The number of new label values 32745480ffeSMatthew G. Knepley - values - The set of new values to relabel keys with 32845480ffeSMatthew G. Knepley 32945480ffeSMatthew G. Knepley Note: This is used internally when boundary label values are specified from the command line. 33045480ffeSMatthew G. Knepley 33145480ffeSMatthew G. Knepley Level: intermediate 33245480ffeSMatthew G. Knepley 333db781477SPatrick Sanan .seealso: `PetscWeakFormReplaceLabel()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 33445480ffeSMatthew G. Knepley @*/ 335d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormRewriteKeys(PetscWeakForm wf, DMLabel label, PetscInt Nv, const PetscInt values[]) 336d71ae5a4SJacob Faibussowitsch { 33706ad1575SMatthew G. Knepley PetscInt f; 33845480ffeSMatthew G. Knepley 33945480ffeSMatthew G. Knepley PetscFunctionBegin; 3409566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormRewriteKeys_Internal(wf, wf->form[f], label, Nv, values)); 34145480ffeSMatthew G. Knepley PetscFunctionReturn(0); 34245480ffeSMatthew G. Knepley } 34345480ffeSMatthew G. Knepley 344d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscWeakFormReplaceLabel_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label) 345d71ae5a4SJacob Faibussowitsch { 346d700741fSMatthew G. Knepley PetscFormKey *keys; 347d700741fSMatthew G. Knepley PetscInt n, i, off = 0, maxFuncs = 0; 348d700741fSMatthew G. Knepley void (**tmpf)(); 349d700741fSMatthew G. Knepley const char *name = NULL; 350d700741fSMatthew G. Knepley 351d700741fSMatthew G. Knepley PetscFunctionBegin; 3529566063dSJacob Faibussowitsch if (label) PetscCall(PetscObjectGetName((PetscObject)label, &name)); 3539566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(hmap, &n)); 3549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &keys)); 3559566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetKeys(hmap, &off, keys)); 356d700741fSMatthew G. Knepley for (i = 0; i < n; ++i) { 357d700741fSMatthew G. Knepley PetscBool match = PETSC_FALSE; 358d700741fSMatthew G. Knepley const char *lname = NULL; 359d700741fSMatthew G. Knepley 360d700741fSMatthew G. Knepley if (label == keys[i].label) continue; 3619566063dSJacob Faibussowitsch if (keys[i].label) PetscCall(PetscObjectGetName((PetscObject)keys[i].label, &lname)); 3629566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(name, lname, &match)); 363d700741fSMatthew G. Knepley if ((!name && !lname) || match) { 364d700741fSMatthew G. Knepley void (**funcs)(); 365d700741fSMatthew G. Knepley PetscInt Nf; 366d700741fSMatthew G. Knepley 3679566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 368d700741fSMatthew G. Knepley maxFuncs = PetscMax(maxFuncs, Nf); 369d700741fSMatthew G. Knepley } 370d700741fSMatthew G. Knepley } 371d700741fSMatthew G. Knepley /* Need temp space because chunk buffer can be reallocated in SetFunction() call */ 3729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxFuncs, &tmpf)); 373d700741fSMatthew G. Knepley for (i = 0; i < n; ++i) { 374d700741fSMatthew G. Knepley PetscBool match = PETSC_FALSE; 375d700741fSMatthew G. Knepley const char *lname = NULL; 376d700741fSMatthew G. Knepley 377d700741fSMatthew G. Knepley if (label == keys[i].label) continue; 3789566063dSJacob Faibussowitsch if (keys[i].label) PetscCall(PetscObjectGetName((PetscObject)keys[i].label, &lname)); 3799566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(name, lname, &match)); 380d700741fSMatthew G. Knepley if ((!name && !lname) || match) { 381d700741fSMatthew G. Knepley void (**funcs)(); 382d700741fSMatthew G. Knepley PetscInt Nf, j; 383d700741fSMatthew G. Knepley 3849566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 385d700741fSMatthew G. Knepley for (j = 0; j < Nf; ++j) tmpf[j] = funcs[j]; 3869566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, label, keys[i].value, keys[i].field, keys[i].part, Nf, tmpf)); 3879566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL)); 388d700741fSMatthew G. Knepley } 389d700741fSMatthew G. Knepley } 3909566063dSJacob Faibussowitsch PetscCall(PetscFree(tmpf)); 3919566063dSJacob Faibussowitsch PetscCall(PetscFree(keys)); 392d700741fSMatthew G. Knepley PetscFunctionReturn(0); 393d700741fSMatthew G. Knepley } 394d700741fSMatthew G. Knepley 395d700741fSMatthew G. Knepley /*@C 396d700741fSMatthew G. Knepley PetscWeakFormReplaceLabel - Change any key on a label of the same name to use the new label 397d700741fSMatthew G. Knepley 398d700741fSMatthew G. Knepley Not Collective 399d700741fSMatthew G. Knepley 400d700741fSMatthew G. Knepley Input Parameters: 401d700741fSMatthew G. Knepley + wf - The original PetscWeakForm 402d700741fSMatthew G. Knepley - label - The label to change keys for 403d700741fSMatthew G. Knepley 404d700741fSMatthew G. Knepley Note: This is used internally when meshes are modified 405d700741fSMatthew G. Knepley 406d700741fSMatthew G. Knepley Level: intermediate 407d700741fSMatthew G. Knepley 408db781477SPatrick Sanan .seealso: `PetscWeakFormRewriteKeys()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 409d700741fSMatthew G. Knepley @*/ 410d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormReplaceLabel(PetscWeakForm wf, DMLabel label) 411d71ae5a4SJacob Faibussowitsch { 412d700741fSMatthew G. Knepley PetscInt f; 413d700741fSMatthew G. Knepley 414d700741fSMatthew G. Knepley PetscFunctionBegin; 4159566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormReplaceLabel_Internal(wf, wf->form[f], label)); 416d700741fSMatthew G. Knepley PetscFunctionReturn(0); 417d700741fSMatthew G. Knepley } 418d700741fSMatthew G. Knepley 419d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormClearIndex(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscWeakFormKind kind, PetscInt ind) 420d71ae5a4SJacob Faibussowitsch { 42106ad1575SMatthew G. Knepley PetscFunctionBegin; 4229566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClearIndexFunction_Private(wf, wf->form[kind], label, val, f, part, ind)); 42306ad1575SMatthew G. Knepley PetscFunctionReturn(0); 42406ad1575SMatthew G. Knepley } 42506ad1575SMatthew G. Knepley 426d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n, void (***obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 427d71ae5a4SJacob Faibussowitsch { 4286528b96dSMatthew G. Knepley PetscFunctionBegin; 4299566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (***)(void))obj)); 4306528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4316528b96dSMatthew G. Knepley } 4326528b96dSMatthew G. Knepley 433d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n, void (**obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 434d71ae5a4SJacob Faibussowitsch { 4356528b96dSMatthew G. Knepley PetscFunctionBegin; 4369566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (**)(void))obj)); 4376528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4386528b96dSMatthew G. Knepley } 4396528b96dSMatthew G. Knepley 440d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormAddObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, void (*obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 441d71ae5a4SJacob Faibussowitsch { 4426528b96dSMatthew G. Knepley PetscFunctionBegin; 4439566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, (void (*)(void))obj)); 4446528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4456528b96dSMatthew G. Knepley } 4466528b96dSMatthew G. Knepley 447d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt ind, void (**obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 448d71ae5a4SJacob Faibussowitsch { 4496528b96dSMatthew G. Knepley PetscFunctionBegin; 4509566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (void (**)(void))obj)); 4516528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4526528b96dSMatthew G. Knepley } 4536528b96dSMatthew G. Knepley 454d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt ind, void (*obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 455d71ae5a4SJacob Faibussowitsch { 4566528b96dSMatthew G. Knepley PetscFunctionBegin; 4579566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (void (*)(void))obj)); 4586528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4596528b96dSMatthew G. Knepley } 4606528b96dSMatthew G. Knepley 461d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n0, void (***f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 462d71ae5a4SJacob Faibussowitsch { 4636528b96dSMatthew G. Knepley PetscFunctionBegin; 4649566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (***)(void))f0)); 4659566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (***)(void))f1)); 4666528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4676528b96dSMatthew G. Knepley } 4686528b96dSMatthew G. Knepley 469d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormAddResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 470d71ae5a4SJacob Faibussowitsch { 4716528b96dSMatthew G. Knepley PetscFunctionBegin; 4729566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, (void (*)(void))f0)); 4739566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, (void (*)(void))f1)); 4746528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4756528b96dSMatthew G. Knepley } 4766528b96dSMatthew G. Knepley 477d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n0, void (**f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 478d71ae5a4SJacob Faibussowitsch { 4796528b96dSMatthew G. Knepley PetscFunctionBegin; 4809566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (**)(void))f0)); 4819566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (**)(void))f1)); 4826528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4836528b96dSMatthew G. Knepley } 4846528b96dSMatthew G. Knepley 485d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetIndexResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt i0, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 486d71ae5a4SJacob Faibussowitsch { 4876528b96dSMatthew G. Knepley PetscFunctionBegin; 4889566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, i0, (void (*)(void))f0)); 4899566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, i1, (void (*)(void))f1)); 4906528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4916528b96dSMatthew G. Knepley } 4926528b96dSMatthew G. Knepley 493d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n0, void (***f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 494d71ae5a4SJacob Faibussowitsch { 4956528b96dSMatthew G. Knepley PetscFunctionBegin; 4969566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (***)(void))f0)); 4979566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (***)(void))f1)); 4986528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4996528b96dSMatthew G. Knepley } 5006528b96dSMatthew G. Knepley 501d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormAddBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 502d71ae5a4SJacob Faibussowitsch { 5036528b96dSMatthew G. Knepley PetscFunctionBegin; 5049566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, (void (*)(void))f0)); 5059566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, (void (*)(void))f1)); 5066528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5076528b96dSMatthew G. Knepley } 5086528b96dSMatthew G. Knepley 509d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n0, void (**f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 510d71ae5a4SJacob Faibussowitsch { 5116528b96dSMatthew G. Knepley PetscFunctionBegin; 5129566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (**)(void))f0)); 5139566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (**)(void))f1)); 5146528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5156528b96dSMatthew G. Knepley } 5166528b96dSMatthew G. Knepley 517d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetIndexBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt i0, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 518d71ae5a4SJacob Faibussowitsch { 5196528b96dSMatthew G. Knepley PetscFunctionBegin; 5209566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, i0, (void (*)(void))f0)); 5219566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, i1, (void (*)(void))f1)); 5226528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5236528b96dSMatthew G. Knepley } 5246528b96dSMatthew G. Knepley 525d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormHasJacobian(PetscWeakForm wf, PetscBool *hasJac) 526d71ae5a4SJacob Faibussowitsch { 5276528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 5286528b96dSMatthew G. Knepley 5296528b96dSMatthew G. Knepley PetscFunctionBegin; 5306528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 5316528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJac, 2); 5329566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G0], &n0)); 5339566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G1], &n1)); 5349566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G2], &n2)); 5359566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G3], &n3)); 5366528b96dSMatthew G. Knepley *hasJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 5376528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5386528b96dSMatthew G. Knepley } 5396528b96dSMatthew G. Knepley 540d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 541d71ae5a4SJacob Faibussowitsch { 5426528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5436528b96dSMatthew G. Knepley 5446528b96dSMatthew G. Knepley PetscFunctionBegin; 5459566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (***)(void))g0)); 5469566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (***)(void))g1)); 5479566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (***)(void))g2)); 5489566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (***)(void))g3)); 5496528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5506528b96dSMatthew G. Knepley } 5516528b96dSMatthew G. Knepley 552d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormAddJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 553d71ae5a4SJacob Faibussowitsch { 5546528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5556528b96dSMatthew G. Knepley 5566528b96dSMatthew G. Knepley PetscFunctionBegin; 5579566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, (void (*)(void))g0)); 5589566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, (void (*)(void))g1)); 5599566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, (void (*)(void))g2)); 5609566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, (void (*)(void))g3)); 5616528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5626528b96dSMatthew G. Knepley } 5636528b96dSMatthew G. Knepley 564d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 565d71ae5a4SJacob Faibussowitsch { 5666528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5676528b96dSMatthew G. Knepley 5686528b96dSMatthew G. Knepley PetscFunctionBegin; 5699566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (**)(void))g0)); 5709566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (**)(void))g1)); 5719566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (**)(void))g2)); 5729566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (**)(void))g3)); 5736528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5746528b96dSMatthew G. Knepley } 5756528b96dSMatthew G. Knepley 576d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetIndexJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 577d71ae5a4SJacob Faibussowitsch { 5786528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5796528b96dSMatthew G. Knepley 5806528b96dSMatthew G. Knepley PetscFunctionBegin; 5819566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, i0, (void (*)(void))g0)); 5829566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, i1, (void (*)(void))g1)); 5839566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, i2, (void (*)(void))g2)); 5849566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, i3, (void (*)(void))g3)); 5856528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5866528b96dSMatthew G. Knepley } 5876528b96dSMatthew G. Knepley 588d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormHasJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre) 589d71ae5a4SJacob Faibussowitsch { 5906528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 5916528b96dSMatthew G. Knepley 5926528b96dSMatthew G. Knepley PetscFunctionBegin; 5936528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 5946528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJacPre, 2); 5959566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP0], &n0)); 5969566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP1], &n1)); 5979566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP2], &n2)); 5989566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP3], &n3)); 5996528b96dSMatthew G. Knepley *hasJacPre = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 6006528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6016528b96dSMatthew G. Knepley } 6026528b96dSMatthew G. Knepley 603d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 604d71ae5a4SJacob Faibussowitsch { 6056528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6066528b96dSMatthew G. Knepley 6076528b96dSMatthew G. Knepley PetscFunctionBegin; 6089566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (***)(void))g0)); 6099566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (***)(void))g1)); 6109566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (***)(void))g2)); 6119566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (***)(void))g3)); 6126528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6136528b96dSMatthew G. Knepley } 6146528b96dSMatthew G. Knepley 615d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormAddJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 616d71ae5a4SJacob Faibussowitsch { 6176528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6186528b96dSMatthew G. Knepley 6196528b96dSMatthew G. Knepley PetscFunctionBegin; 6209566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, (void (*)(void))g0)); 6219566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, (void (*)(void))g1)); 6229566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, (void (*)(void))g2)); 6239566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, (void (*)(void))g3)); 6246528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6256528b96dSMatthew G. Knepley } 6266528b96dSMatthew G. Knepley 627d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 628d71ae5a4SJacob Faibussowitsch { 6296528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6306528b96dSMatthew G. Knepley 6316528b96dSMatthew G. Knepley PetscFunctionBegin; 6329566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (**)(void))g0)); 6339566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (**)(void))g1)); 6349566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (**)(void))g2)); 6359566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (**)(void))g3)); 6366528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6376528b96dSMatthew G. Knepley } 6386528b96dSMatthew G. Knepley 639d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetIndexJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 640d71ae5a4SJacob Faibussowitsch { 6416528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6426528b96dSMatthew G. Knepley 6436528b96dSMatthew G. Knepley PetscFunctionBegin; 6449566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, i0, (void (*)(void))g0)); 6459566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, i1, (void (*)(void))g1)); 6469566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, i2, (void (*)(void))g2)); 6479566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, i3, (void (*)(void))g3)); 6486528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6496528b96dSMatthew G. Knepley } 6506528b96dSMatthew G. Knepley 651d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormHasBdJacobian(PetscWeakForm wf, PetscBool *hasJac) 652d71ae5a4SJacob Faibussowitsch { 6536528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 6546528b96dSMatthew G. Knepley 6556528b96dSMatthew G. Knepley PetscFunctionBegin; 6566528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 6576528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJac, 2); 6589566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG0], &n0)); 6599566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG1], &n1)); 6609566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG2], &n2)); 6619566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG3], &n3)); 6626528b96dSMatthew G. Knepley *hasJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 6636528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6646528b96dSMatthew G. Knepley } 6656528b96dSMatthew G. Knepley 666d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 667d71ae5a4SJacob Faibussowitsch { 6686528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6696528b96dSMatthew G. Knepley 6706528b96dSMatthew G. Knepley PetscFunctionBegin; 6719566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (***)(void))g0)); 6729566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (***)(void))g1)); 6739566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (***)(void))g2)); 6749566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (***)(void))g3)); 6756528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6766528b96dSMatthew G. Knepley } 6776528b96dSMatthew G. Knepley 678d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormAddBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 679d71ae5a4SJacob Faibussowitsch { 6806528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6816528b96dSMatthew G. Knepley 6826528b96dSMatthew G. Knepley PetscFunctionBegin; 6839566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, (void (*)(void))g0)); 6849566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, (void (*)(void))g1)); 6859566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, (void (*)(void))g2)); 6869566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, (void (*)(void))g3)); 6876528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6886528b96dSMatthew G. Knepley } 6896528b96dSMatthew G. Knepley 690d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 691d71ae5a4SJacob Faibussowitsch { 6926528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6936528b96dSMatthew G. Knepley 6946528b96dSMatthew G. Knepley PetscFunctionBegin; 6959566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (**)(void))g0)); 6969566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (**)(void))g1)); 6979566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (**)(void))g2)); 6989566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (**)(void))g3)); 6996528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7006528b96dSMatthew G. Knepley } 7016528b96dSMatthew G. Knepley 702d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetIndexBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 703d71ae5a4SJacob Faibussowitsch { 7046528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7056528b96dSMatthew G. Knepley 7066528b96dSMatthew G. Knepley PetscFunctionBegin; 7079566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, i0, (void (*)(void))g0)); 7089566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, i1, (void (*)(void))g1)); 7099566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, i2, (void (*)(void))g2)); 7109566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, i3, (void (*)(void))g3)); 7116528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7126528b96dSMatthew G. Knepley } 7136528b96dSMatthew G. Knepley 714d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre) 715d71ae5a4SJacob Faibussowitsch { 7166528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 7176528b96dSMatthew G. Knepley 7186528b96dSMatthew G. Knepley PetscFunctionBegin; 7196528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 7206528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJacPre, 2); 7219566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP0], &n0)); 7229566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP1], &n1)); 7239566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP2], &n2)); 7249566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP3], &n3)); 7256528b96dSMatthew G. Knepley *hasJacPre = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 7266528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7276528b96dSMatthew G. Knepley } 7286528b96dSMatthew G. Knepley 729d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 730d71ae5a4SJacob Faibussowitsch { 7316528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7326528b96dSMatthew G. Knepley 7336528b96dSMatthew G. Knepley PetscFunctionBegin; 7349566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (***)(void))g0)); 7359566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (***)(void))g1)); 7369566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (***)(void))g2)); 7379566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (***)(void))g3)); 7386528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7396528b96dSMatthew G. Knepley } 7406528b96dSMatthew G. Knepley 741d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormAddBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 742d71ae5a4SJacob Faibussowitsch { 7436528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7446528b96dSMatthew G. Knepley 7456528b96dSMatthew G. Knepley PetscFunctionBegin; 7469566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, (void (*)(void))g0)); 7479566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, (void (*)(void))g1)); 7489566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, (void (*)(void))g2)); 7499566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, (void (*)(void))g3)); 7506528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7516528b96dSMatthew G. Knepley } 7526528b96dSMatthew G. Knepley 753d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 754d71ae5a4SJacob Faibussowitsch { 7556528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7566528b96dSMatthew G. Knepley 7576528b96dSMatthew G. Knepley PetscFunctionBegin; 7589566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (**)(void))g0)); 7599566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (**)(void))g1)); 7609566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (**)(void))g2)); 7619566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (**)(void))g3)); 7626528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7636528b96dSMatthew G. Knepley } 7646528b96dSMatthew G. Knepley 765d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetIndexBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 766d71ae5a4SJacob Faibussowitsch { 7676528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7686528b96dSMatthew G. Knepley 7696528b96dSMatthew G. Knepley PetscFunctionBegin; 7709566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, i0, (void (*)(void))g0)); 7719566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, i1, (void (*)(void))g1)); 7729566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, i2, (void (*)(void))g2)); 7739566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, i3, (void (*)(void))g3)); 7746528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7756528b96dSMatthew G. Knepley } 7766528b96dSMatthew G. Knepley 777d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormHasDynamicJacobian(PetscWeakForm wf, PetscBool *hasDynJac) 778d71ae5a4SJacob Faibussowitsch { 7796528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 7806528b96dSMatthew G. Knepley 7816528b96dSMatthew G. Knepley PetscFunctionBegin; 7826528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 7836528b96dSMatthew G. Knepley PetscValidBoolPointer(hasDynJac, 2); 7849566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT0], &n0)); 7859566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT1], &n1)); 7869566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT2], &n2)); 7879566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT3], &n3)); 7886528b96dSMatthew G. Knepley *hasDynJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 7896528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7906528b96dSMatthew G. Knepley } 7916528b96dSMatthew G. Knepley 792d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 793d71ae5a4SJacob Faibussowitsch { 7946528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7956528b96dSMatthew G. Knepley 7966528b96dSMatthew G. Knepley PetscFunctionBegin; 7979566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (***)(void))g0)); 7989566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (***)(void))g1)); 7999566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (***)(void))g2)); 8009566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (***)(void))g3)); 8016528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8026528b96dSMatthew G. Knepley } 8036528b96dSMatthew G. Knepley 804d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormAddDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 805d71ae5a4SJacob Faibussowitsch { 8066528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 8076528b96dSMatthew G. Knepley 8086528b96dSMatthew G. Knepley PetscFunctionBegin; 8099566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, (void (*)(void))g0)); 8109566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, (void (*)(void))g1)); 8119566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, (void (*)(void))g2)); 8129566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, (void (*)(void))g3)); 8136528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8146528b96dSMatthew G. Knepley } 8156528b96dSMatthew G. Knepley 816d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 817d71ae5a4SJacob Faibussowitsch { 8186528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 8196528b96dSMatthew G. Knepley 8206528b96dSMatthew G. Knepley PetscFunctionBegin; 8219566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (**)(void))g0)); 8229566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (**)(void))g1)); 8239566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (**)(void))g2)); 8249566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (**)(void))g3)); 8256528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8266528b96dSMatthew G. Knepley } 8276528b96dSMatthew G. Knepley 828d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetIndexDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) 829d71ae5a4SJacob Faibussowitsch { 8306528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 8316528b96dSMatthew G. Knepley 8326528b96dSMatthew G. Knepley PetscFunctionBegin; 8339566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, i0, (void (*)(void))g0)); 8349566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, i1, (void (*)(void))g1)); 8359566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, i2, (void (*)(void))g2)); 8369566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, i3, (void (*)(void))g3)); 8376528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8386528b96dSMatthew G. Knepley } 8396528b96dSMatthew G. Knepley 840d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n, void (***r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)) 841d71ae5a4SJacob Faibussowitsch { 8426528b96dSMatthew G. Knepley PetscFunctionBegin; 8439566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (***)(void))r)); 8446528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8456528b96dSMatthew G. Knepley } 8466528b96dSMatthew G. Knepley 847d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n, void (**r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)) 848d71ae5a4SJacob Faibussowitsch { 8496528b96dSMatthew G. Knepley PetscFunctionBegin; 8509566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (**)(void))r)); 8516528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8526528b96dSMatthew G. Knepley } 8536528b96dSMatthew G. Knepley 854d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetIndexRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt i, void (*r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)) 855d71ae5a4SJacob Faibussowitsch { 8566528b96dSMatthew G. Knepley PetscFunctionBegin; 8579566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, i, (void (*)(void))r)); 8586528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8596528b96dSMatthew G. Knepley } 8606528b96dSMatthew G. Knepley 8616528b96dSMatthew G. Knepley /*@ 8626528b96dSMatthew G. Knepley PetscWeakFormGetNumFields - Returns the number of fields 8636528b96dSMatthew G. Knepley 8646528b96dSMatthew G. Knepley Not collective 8656528b96dSMatthew G. Knepley 8666528b96dSMatthew G. Knepley Input Parameter: 8676528b96dSMatthew G. Knepley . wf - The PetscWeakForm object 8686528b96dSMatthew G. Knepley 8696528b96dSMatthew G. Knepley Output Parameter: 870a5b23f4aSJose E. Roman . Nf - The number of fields 8716528b96dSMatthew G. Knepley 8726528b96dSMatthew G. Knepley Level: beginner 8736528b96dSMatthew G. Knepley 874db781477SPatrick Sanan .seealso: `PetscWeakFormSetNumFields()`, `PetscWeakFormCreate()` 8756528b96dSMatthew G. Knepley @*/ 876d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormGetNumFields(PetscWeakForm wf, PetscInt *Nf) 877d71ae5a4SJacob Faibussowitsch { 8786528b96dSMatthew G. Knepley PetscFunctionBegin; 8796528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 880dadcf809SJacob Faibussowitsch PetscValidIntPointer(Nf, 2); 8816528b96dSMatthew G. Knepley *Nf = wf->Nf; 8826528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8836528b96dSMatthew G. Knepley } 8846528b96dSMatthew G. Knepley 8856528b96dSMatthew G. Knepley /*@ 8866528b96dSMatthew G. Knepley PetscWeakFormSetNumFields - Sets the number of fields 8876528b96dSMatthew G. Knepley 8886528b96dSMatthew G. Knepley Not collective 8896528b96dSMatthew G. Knepley 8906528b96dSMatthew G. Knepley Input Parameters: 8916528b96dSMatthew G. Knepley + wf - The PetscWeakForm object 8926528b96dSMatthew G. Knepley - Nf - The number of fields 8936528b96dSMatthew G. Knepley 8946528b96dSMatthew G. Knepley Level: beginner 8956528b96dSMatthew G. Knepley 896db781477SPatrick Sanan .seealso: `PetscWeakFormGetNumFields()`, `PetscWeakFormCreate()` 8976528b96dSMatthew G. Knepley @*/ 898d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormSetNumFields(PetscWeakForm wf, PetscInt Nf) 899d71ae5a4SJacob Faibussowitsch { 9006528b96dSMatthew G. Knepley PetscFunctionBegin; 9016528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 9026528b96dSMatthew G. Knepley wf->Nf = Nf; 9036528b96dSMatthew G. Knepley PetscFunctionReturn(0); 9046528b96dSMatthew G. Knepley } 9056528b96dSMatthew G. Knepley 9066528b96dSMatthew G. Knepley /*@ 9076528b96dSMatthew G. Knepley PetscWeakFormDestroy - Destroys a PetscWeakForm object 9086528b96dSMatthew G. Knepley 9096528b96dSMatthew G. Knepley Collective on wf 9106528b96dSMatthew G. Knepley 9116528b96dSMatthew G. Knepley Input Parameter: 9126528b96dSMatthew G. Knepley . wf - the PetscWeakForm object to destroy 9136528b96dSMatthew G. Knepley 9146528b96dSMatthew G. Knepley Level: developer 9156528b96dSMatthew G. Knepley 916db781477SPatrick Sanan .seealso `PetscWeakFormCreate()`, `PetscWeakFormView()` 9176528b96dSMatthew G. Knepley @*/ 918d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormDestroy(PetscWeakForm *wf) 919d71ae5a4SJacob Faibussowitsch { 92006ad1575SMatthew G. Knepley PetscInt f; 9216528b96dSMatthew G. Knepley 9226528b96dSMatthew G. Knepley PetscFunctionBegin; 9236528b96dSMatthew G. Knepley if (!*wf) PetscFunctionReturn(0); 9246528b96dSMatthew G. Knepley PetscValidHeaderSpecific((*wf), PETSCWEAKFORM_CLASSID, 1); 9256528b96dSMatthew G. Knepley 9269371c9d4SSatish Balay if (--((PetscObject)(*wf))->refct > 0) { 9279371c9d4SSatish Balay *wf = NULL; 9289371c9d4SSatish Balay PetscFunctionReturn(0); 9299371c9d4SSatish Balay } 9306528b96dSMatthew G. Knepley ((PetscObject)(*wf))->refct = 0; 9319566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferDestroy(&(*wf)->funcs)); 9329566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormDestroy(&(*wf)->form[f])); 9339566063dSJacob Faibussowitsch PetscCall(PetscFree((*wf)->form)); 9349566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(wf)); 9356528b96dSMatthew G. Knepley PetscFunctionReturn(0); 9366528b96dSMatthew G. Knepley } 9376528b96dSMatthew G. Knepley 938d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscWeakFormViewTable_Ascii(PetscWeakForm wf, PetscViewer viewer, PetscBool splitField, const char tableName[], PetscHMapForm map) 939d71ae5a4SJacob Faibussowitsch { 94045480ffeSMatthew G. Knepley PetscInt Nf = wf->Nf, Nk, k; 9416528b96dSMatthew G. Knepley 9426528b96dSMatthew G. Knepley PetscFunctionBegin; 9439566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(map, &Nk)); 9446528b96dSMatthew G. Knepley if (Nk) { 94506ad1575SMatthew G. Knepley PetscFormKey *keys; 946165f9cc3SJed Brown void (**funcs)(void) = NULL; 9475fedec97SMatthew G. Knepley const char **names; 9485fedec97SMatthew G. Knepley PetscInt *values, *idx1, *idx2, *idx; 9495fedec97SMatthew G. Knepley PetscBool showPart = PETSC_FALSE, showPointer = PETSC_FALSE; 9505fedec97SMatthew G. Knepley PetscInt off = 0; 9516528b96dSMatthew G. Knepley 9529566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(Nk, &keys, Nk, &names, Nk, &values, Nk, &idx1, Nk, &idx2, Nk, &idx)); 9539566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetKeys(map, &off, keys)); 9545fedec97SMatthew G. Knepley /* Sort keys by label name and value */ 9555fedec97SMatthew G. Knepley { 9565fedec97SMatthew G. Knepley /* First sort values */ 9579371c9d4SSatish Balay for (k = 0; k < Nk; ++k) { 9589371c9d4SSatish Balay values[k] = keys[k].value; 9599371c9d4SSatish Balay idx1[k] = k; 9609371c9d4SSatish Balay } 9619566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(Nk, values, idx1)); 9625fedec97SMatthew G. Knepley /* If the string sort is stable, it will be sorted correctly overall */ 9635fedec97SMatthew G. Knepley for (k = 0; k < Nk; ++k) { 9649566063dSJacob Faibussowitsch if (keys[idx1[k]].label) PetscCall(PetscObjectGetName((PetscObject)keys[idx1[k]].label, &names[k])); 965ad540459SPierre Jolivet else names[k] = ""; 9665fedec97SMatthew G. Knepley idx2[k] = k; 9675fedec97SMatthew G. Knepley } 9689566063dSJacob Faibussowitsch PetscCall(PetscSortStrWithPermutation(Nk, names, idx2)); 9695fedec97SMatthew G. Knepley for (k = 0; k < Nk; ++k) { 9709566063dSJacob Faibussowitsch if (keys[k].label) PetscCall(PetscObjectGetName((PetscObject)keys[k].label, &names[k])); 971ad540459SPierre Jolivet else names[k] = ""; 9725fedec97SMatthew G. Knepley idx[k] = idx1[idx2[k]]; 9735fedec97SMatthew G. Knepley } 9745fedec97SMatthew G. Knepley } 9759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", tableName)); 9769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 9776528b96dSMatthew G. Knepley for (k = 0; k < Nk; ++k) { 9780944a4d6SMatthew G. Knepley if (keys[k].part != 0) showPart = PETSC_TRUE; 9790944a4d6SMatthew G. Knepley } 9800944a4d6SMatthew G. Knepley for (k = 0; k < Nk; ++k) { 9815fedec97SMatthew G. Knepley const PetscInt i = idx[k]; 9825fedec97SMatthew G. Knepley PetscInt n, f; 9835fedec97SMatthew G. Knepley 9845fedec97SMatthew G. Knepley if (keys[i].label) { 985*c75bfeddSPierre Jolivet if (showPointer) PetscCall(PetscViewerASCIIPrintf(viewer, "(%s:%p, %" PetscInt_FMT ") ", names[i], (void *)keys[i].label, keys[i].value)); 98663a3b9bcSJacob Faibussowitsch else PetscCall(PetscViewerASCIIPrintf(viewer, "(%s, %" PetscInt_FMT ") ", names[i], keys[i].value)); 98763a3b9bcSJacob Faibussowitsch } 9889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 98963a3b9bcSJacob Faibussowitsch if (splitField) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ", %" PetscInt_FMT ") ", keys[i].field / Nf, keys[i].field % Nf)); 99063a3b9bcSJacob Faibussowitsch else PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") ", keys[i].field)); 99163a3b9bcSJacob Faibussowitsch if (showPart) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") ", keys[i].part)); 9929566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, map, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &n, &funcs)); 9935fedec97SMatthew G. Knepley for (f = 0; f < n; ++f) { 994258ec3d2SMatthew G. Knepley char *fname; 9955fedec97SMatthew G. Knepley size_t len, l; 996258ec3d2SMatthew G. Knepley 9979566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 9989566063dSJacob Faibussowitsch PetscCall(PetscDLAddr(funcs[f], &fname)); 9995fedec97SMatthew G. Knepley if (fname) { 10005fedec97SMatthew G. Knepley /* Eliminate argument types */ 10019566063dSJacob Faibussowitsch PetscCall(PetscStrlen(fname, &len)); 10029371c9d4SSatish Balay for (l = 0; l < len; ++l) 10039371c9d4SSatish Balay if (fname[l] == '(') { 10049371c9d4SSatish Balay fname[l] = '\0'; 10059371c9d4SSatish Balay break; 10069371c9d4SSatish Balay } 10079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%s", fname)); 10085fedec97SMatthew G. Knepley } else if (showPointer) { 10099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%p", funcs[f])); 10105fedec97SMatthew G. Knepley } 10119566063dSJacob Faibussowitsch PetscCall(PetscFree(fname)); 10126528b96dSMatthew G. Knepley } 10139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 10149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 10156528b96dSMatthew G. Knepley } 10169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 10179566063dSJacob Faibussowitsch PetscCall(PetscFree6(keys, names, values, idx1, idx2, idx)); 10186528b96dSMatthew G. Knepley } 10196528b96dSMatthew G. Knepley PetscFunctionReturn(0); 10206528b96dSMatthew G. Knepley } 10216528b96dSMatthew G. Knepley 1022d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscWeakFormView_Ascii(PetscWeakForm wf, PetscViewer viewer) 1023d71ae5a4SJacob Faibussowitsch { 10246528b96dSMatthew G. Knepley PetscViewerFormat format; 102506ad1575SMatthew G. Knepley PetscInt f; 10266528b96dSMatthew G. Knepley 10276528b96dSMatthew G. Knepley PetscFunctionBegin; 10289566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 102963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Weak Form System with %" PetscInt_FMT " fields\n", wf->Nf)); 10309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 103148a46eb9SPierre Jolivet for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormViewTable_Ascii(wf, viewer, PETSC_TRUE, PetscWeakFormKinds[f], wf->form[f])); 10329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 10336528b96dSMatthew G. Knepley PetscFunctionReturn(0); 10346528b96dSMatthew G. Knepley } 10356528b96dSMatthew G. Knepley 10366528b96dSMatthew G. Knepley /*@C 10376528b96dSMatthew G. Knepley PetscWeakFormView - Views a PetscWeakForm 10386528b96dSMatthew G. Knepley 10396528b96dSMatthew G. Knepley Collective on wf 10406528b96dSMatthew G. Knepley 1041d8d19677SJose E. Roman Input Parameters: 10426528b96dSMatthew G. Knepley + wf - the PetscWeakForm object to view 10436528b96dSMatthew G. Knepley - v - the viewer 10446528b96dSMatthew G. Knepley 10456528b96dSMatthew G. Knepley Level: developer 10466528b96dSMatthew G. Knepley 1047db781477SPatrick Sanan .seealso `PetscWeakFormDestroy()`, `PetscWeakFormCreate()` 10486528b96dSMatthew G. Knepley @*/ 1049d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormView(PetscWeakForm wf, PetscViewer v) 1050d71ae5a4SJacob Faibussowitsch { 10516528b96dSMatthew G. Knepley PetscBool iascii; 10526528b96dSMatthew G. Knepley 10536528b96dSMatthew G. Knepley PetscFunctionBegin; 10546528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 10559566063dSJacob Faibussowitsch if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)wf), &v)); 1056ad540459SPierre Jolivet else PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2); 10579566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii)); 10589566063dSJacob Faibussowitsch if (iascii) PetscCall(PetscWeakFormView_Ascii(wf, v)); 1059dbbe0bcdSBarry Smith PetscTryTypeMethod(wf, view, v); 10606528b96dSMatthew G. Knepley PetscFunctionReturn(0); 10616528b96dSMatthew G. Knepley } 10626528b96dSMatthew G. Knepley 10636528b96dSMatthew G. Knepley /*@ 10646528b96dSMatthew G. Knepley PetscWeakFormCreate - Creates an empty PetscWeakForm object. 10656528b96dSMatthew G. Knepley 10666528b96dSMatthew G. Knepley Collective 10676528b96dSMatthew G. Knepley 10686528b96dSMatthew G. Knepley Input Parameter: 10696528b96dSMatthew G. Knepley . comm - The communicator for the PetscWeakForm object 10706528b96dSMatthew G. Knepley 10716528b96dSMatthew G. Knepley Output Parameter: 10726528b96dSMatthew G. Knepley . wf - The PetscWeakForm object 10736528b96dSMatthew G. Knepley 10746528b96dSMatthew G. Knepley Level: beginner 10756528b96dSMatthew G. Knepley 1076db781477SPatrick Sanan .seealso: `PetscDS`, `PetscWeakFormDestroy()` 10776528b96dSMatthew G. Knepley @*/ 1078d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscWeakFormCreate(MPI_Comm comm, PetscWeakForm *wf) 1079d71ae5a4SJacob Faibussowitsch { 10806528b96dSMatthew G. Knepley PetscWeakForm p; 108106ad1575SMatthew G. Knepley PetscInt f; 10826528b96dSMatthew G. Knepley 10836528b96dSMatthew G. Knepley PetscFunctionBegin; 10846528b96dSMatthew G. Knepley PetscValidPointer(wf, 2); 10856528b96dSMatthew G. Knepley *wf = NULL; 10869566063dSJacob Faibussowitsch PetscCall(PetscDSInitializePackage()); 10876528b96dSMatthew G. Knepley 10889566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(p, PETSCWEAKFORM_CLASSID, "PetscWeakForm", "Weak Form System", "PetscWeakForm", comm, PetscWeakFormDestroy, PetscWeakFormView)); 10896528b96dSMatthew G. Knepley 10906528b96dSMatthew G. Knepley p->Nf = 0; 10919566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreate(sizeof(&PetscWeakFormCreate), 2, &p->funcs)); 10929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PETSC_NUM_WF, &p->form)); 10939566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormCreate(&p->form[f])); 10946528b96dSMatthew G. Knepley *wf = p; 10956528b96dSMatthew G. Knepley PetscFunctionReturn(0); 10966528b96dSMatthew G. Knepley } 1097