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 79371c9d4SSatish Balay static PetscErrorCode PetscChunkBufferCreate(size_t unitbytes, size_t expected, PetscChunkBuffer **buffer) { 86528b96dSMatthew G. Knepley PetscFunctionBegin; 99566063dSJacob Faibussowitsch PetscCall(PetscNew(buffer)); 109566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(expected * unitbytes, &(*buffer)->array)); 116528b96dSMatthew G. Knepley (*buffer)->size = expected; 126528b96dSMatthew G. Knepley (*buffer)->unitbytes = unitbytes; 136528b96dSMatthew G. Knepley (*buffer)->alloc = expected * unitbytes; 146528b96dSMatthew G. Knepley PetscFunctionReturn(0); 156528b96dSMatthew G. Knepley } 166528b96dSMatthew G. Knepley 179371c9d4SSatish Balay static PetscErrorCode PetscChunkBufferDuplicate(PetscChunkBuffer *buffer, PetscChunkBuffer **bufferNew) { 1845480ffeSMatthew G. Knepley PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall(PetscNew(bufferNew)); 209566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(buffer->size * buffer->unitbytes, &(*bufferNew)->array)); 219566063dSJacob Faibussowitsch PetscCall(PetscMemcpy((*bufferNew)->array, buffer->array, buffer->size * buffer->unitbytes)); 2245480ffeSMatthew G. Knepley (*bufferNew)->size = buffer->size; 2345480ffeSMatthew G. Knepley (*bufferNew)->unitbytes = buffer->unitbytes; 2445480ffeSMatthew G. Knepley (*bufferNew)->alloc = buffer->size * buffer->unitbytes; 2545480ffeSMatthew G. Knepley PetscFunctionReturn(0); 2645480ffeSMatthew G. Knepley } 2745480ffeSMatthew G. Knepley 289371c9d4SSatish Balay static PetscErrorCode PetscChunkBufferDestroy(PetscChunkBuffer **buffer) { 296528b96dSMatthew G. Knepley PetscFunctionBegin; 309566063dSJacob Faibussowitsch PetscCall(PetscFree((*buffer)->array)); 319566063dSJacob Faibussowitsch PetscCall(PetscFree(*buffer)); 326528b96dSMatthew G. Knepley PetscFunctionReturn(0); 336528b96dSMatthew G. Knepley } 346528b96dSMatthew G. Knepley 359371c9d4SSatish Balay static PetscErrorCode PetscChunkBufferCreateChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk) { 366528b96dSMatthew G. Knepley PetscFunctionBegin; 376528b96dSMatthew G. Knepley if ((buffer->size + size) * buffer->unitbytes > buffer->alloc) { 386528b96dSMatthew G. Knepley char *tmp; 396528b96dSMatthew G. Knepley 406528b96dSMatthew G. Knepley if (!buffer->alloc) buffer->alloc = (buffer->size + size) * buffer->unitbytes; 416528b96dSMatthew G. Knepley while ((buffer->size + size) * buffer->unitbytes > buffer->alloc) buffer->alloc *= 2; 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc(buffer->alloc, &tmp)); 439566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(tmp, buffer->array, buffer->size * buffer->unitbytes)); 449566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer->array)); 456528b96dSMatthew G. Knepley buffer->array = tmp; 466528b96dSMatthew G. Knepley } 472baeeceeSMatthew G. Knepley chunk->start = buffer->size * buffer->unitbytes; 486528b96dSMatthew G. Knepley chunk->size = size; 496528b96dSMatthew G. Knepley chunk->reserved = size; 506528b96dSMatthew G. Knepley buffer->size += size; 516528b96dSMatthew G. Knepley PetscFunctionReturn(0); 526528b96dSMatthew G. Knepley } 536528b96dSMatthew G. Knepley 549371c9d4SSatish Balay static PetscErrorCode PetscChunkBufferEnlargeChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk) { 556528b96dSMatthew G. Knepley size_t siz = size; 566528b96dSMatthew G. Knepley 576528b96dSMatthew G. Knepley PetscFunctionBegin; 586528b96dSMatthew G. Knepley if (chunk->size + size > chunk->reserved) { 596528b96dSMatthew G. Knepley PetscChunk newchunk; 606528b96dSMatthew G. Knepley PetscInt reserved = chunk->size; 616528b96dSMatthew G. Knepley 626528b96dSMatthew G. Knepley /* TODO Here if we had a chunk list, we could update them all to reclaim unused space */ 636528b96dSMatthew G. Knepley while (reserved < chunk->size + size) reserved *= 2; 649566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(buffer, (size_t)reserved, &newchunk)); 656528b96dSMatthew G. Knepley newchunk.size = chunk->size + size; 669566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(&buffer->array[newchunk.start], &buffer->array[chunk->start], chunk->size * buffer->unitbytes)); 676528b96dSMatthew G. Knepley *chunk = newchunk; 686528b96dSMatthew G. Knepley } else { 696528b96dSMatthew G. Knepley chunk->size += siz; 706528b96dSMatthew G. Knepley } 716528b96dSMatthew G. Knepley PetscFunctionReturn(0); 726528b96dSMatthew G. Knepley } 736528b96dSMatthew G. Knepley 746528b96dSMatthew G. Knepley /*@C 7506ad1575SMatthew G. Knepley PetscFormKeySort - Sorts an array of PetscFormKey in place in increasing order. 766528b96dSMatthew G. Knepley 776528b96dSMatthew G. Knepley Not Collective 786528b96dSMatthew G. Knepley 796528b96dSMatthew G. Knepley Input Parameters: 806528b96dSMatthew G. Knepley + n - number of values 8106ad1575SMatthew G. Knepley - X - array of PetscFormKey 826528b96dSMatthew G. Knepley 836528b96dSMatthew G. Knepley Level: intermediate 846528b96dSMatthew G. Knepley 85db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()` 866528b96dSMatthew G. Knepley @*/ 879371c9d4SSatish Balay PetscErrorCode PetscFormKeySort(PetscInt n, PetscFormKey arr[]) { 886528b96dSMatthew G. Knepley PetscFunctionBegin; 896528b96dSMatthew G. Knepley if (n <= 1) PetscFunctionReturn(0); 906528b96dSMatthew G. Knepley PetscValidPointer(arr, 2); 919566063dSJacob Faibussowitsch PetscCall(PetscTimSort(n, arr, sizeof(PetscFormKey), Compare_PetscFormKey_Private, NULL)); 926528b96dSMatthew G. Knepley PetscFunctionReturn(0); 936528b96dSMatthew G. Knepley } 946528b96dSMatthew G. Knepley 959371c9d4SSatish Balay PetscErrorCode PetscWeakFormGetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt *n, void (***func)()) { 9606ad1575SMatthew G. Knepley PetscFormKey key; 976528b96dSMatthew G. Knepley PetscChunk chunk; 986528b96dSMatthew G. Knepley 996528b96dSMatthew G. Knepley PetscFunctionBegin; 1009371c9d4SSatish Balay key.label = label; 1019371c9d4SSatish Balay key.value = value; 1029371c9d4SSatish Balay key.field = f; 1039371c9d4SSatish Balay key.part = part; 1049566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1059371c9d4SSatish Balay if (chunk.size < 0) { 1069371c9d4SSatish Balay *n = 0; 1079371c9d4SSatish Balay *func = NULL; 1089371c9d4SSatish Balay } else { 1099371c9d4SSatish Balay *n = chunk.size; 1109371c9d4SSatish Balay *func = (void (**)()) & wf->funcs->array[chunk.start]; 1119371c9d4SSatish Balay } 1126528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1136528b96dSMatthew G. Knepley } 1146528b96dSMatthew G. Knepley 1156528b96dSMatthew G. Knepley /* A NULL argument for func causes this to clear the key */ 1169371c9d4SSatish Balay PetscErrorCode PetscWeakFormSetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt n, void (**func)()) { 11706ad1575SMatthew G. Knepley PetscFormKey key; 1186528b96dSMatthew G. Knepley PetscChunk chunk; 1196528b96dSMatthew G. Knepley PetscInt i; 1206528b96dSMatthew G. Knepley 1216528b96dSMatthew G. Knepley PetscFunctionBegin; 1229371c9d4SSatish Balay key.label = label; 1239371c9d4SSatish Balay key.value = value; 1249371c9d4SSatish Balay key.field = f; 1259371c9d4SSatish Balay key.part = part; 1266528b96dSMatthew G. Knepley if (!func) { 1279566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDel(ht, key)); 1286528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1291baa6e33SBarry Smith } else PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1306528b96dSMatthew G. Knepley if (chunk.size < 0) { 1319566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(wf->funcs, n, &chunk)); 1329566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1336528b96dSMatthew G. Knepley } else if (chunk.size <= n) { 1349566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, n - chunk.size, &chunk)); 1359566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1366528b96dSMatthew G. Knepley } 1372baeeceeSMatthew G. Knepley for (i = 0; i < n; ++i) ((void (**)()) & wf->funcs->array[chunk.start])[i] = func[i]; 1386528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1396528b96dSMatthew G. Knepley } 1406528b96dSMatthew G. Knepley 1419371c9d4SSatish Balay PetscErrorCode PetscWeakFormAddFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, void (*func)()) { 14206ad1575SMatthew G. Knepley PetscFormKey key; 1436528b96dSMatthew G. Knepley PetscChunk chunk; 1446528b96dSMatthew G. Knepley 1456528b96dSMatthew G. Knepley PetscFunctionBegin; 1466528b96dSMatthew G. Knepley if (!func) PetscFunctionReturn(0); 1479371c9d4SSatish Balay key.label = label; 1489371c9d4SSatish Balay key.value = value; 1499371c9d4SSatish Balay key.field = f; 1509371c9d4SSatish Balay key.part = part; 1519566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1526528b96dSMatthew G. Knepley if (chunk.size < 0) { 1539566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(wf->funcs, 1, &chunk)); 1549566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1552baeeceeSMatthew G. Knepley ((void (**)()) & wf->funcs->array[chunk.start])[0] = func; 1566528b96dSMatthew G. Knepley } else { 1579566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, 1, &chunk)); 1589566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1592baeeceeSMatthew G. Knepley ((void (**)()) & wf->funcs->array[chunk.start])[chunk.size - 1] = func; 1606528b96dSMatthew G. Knepley } 1616528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1626528b96dSMatthew G. Knepley } 1636528b96dSMatthew G. Knepley 1649371c9d4SSatish Balay PetscErrorCode PetscWeakFormGetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, void (**func)()) { 16506ad1575SMatthew G. Knepley PetscFormKey key; 1666528b96dSMatthew G. Knepley PetscChunk chunk; 1676528b96dSMatthew G. Knepley 1686528b96dSMatthew G. Knepley PetscFunctionBegin; 1699371c9d4SSatish Balay key.label = label; 1709371c9d4SSatish Balay key.value = value; 1719371c9d4SSatish Balay key.field = f; 1729371c9d4SSatish Balay key.part = part; 1739566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1749371c9d4SSatish Balay if (chunk.size < 0) { 1759371c9d4SSatish Balay *func = NULL; 1769371c9d4SSatish Balay } else { 17763a3b9bcSJacob Faibussowitsch PetscCheck(ind < chunk.size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", ind, chunk.size); 1782baeeceeSMatthew G. Knepley *func = ((void (**)()) & wf->funcs->array[chunk.start])[ind]; 1796528b96dSMatthew G. Knepley } 1806528b96dSMatthew G. Knepley PetscFunctionReturn(0); 1816528b96dSMatthew G. Knepley } 1826528b96dSMatthew G. Knepley 18306ad1575SMatthew G. Knepley /* Ignore a NULL func */ 1849371c9d4SSatish Balay PetscErrorCode PetscWeakFormSetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, void (*func)()) { 18506ad1575SMatthew G. Knepley PetscFormKey key; 1866528b96dSMatthew G. Knepley PetscChunk chunk; 1876528b96dSMatthew G. Knepley 1886528b96dSMatthew G. Knepley PetscFunctionBegin; 18906ad1575SMatthew G. Knepley if (!func) PetscFunctionReturn(0); 1909371c9d4SSatish Balay key.label = label; 1919371c9d4SSatish Balay key.value = value; 1929371c9d4SSatish Balay key.field = f; 1939371c9d4SSatish Balay key.part = part; 1949566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 1956528b96dSMatthew G. Knepley if (chunk.size < 0) { 1969566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreateChunk(wf->funcs, ind + 1, &chunk)); 1979566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 1986528b96dSMatthew G. Knepley } else if (chunk.size <= ind) { 1999566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, ind - chunk.size + 1, &chunk)); 2009566063dSJacob Faibussowitsch PetscCall(PetscHMapFormSet(ht, key, chunk)); 2016528b96dSMatthew G. Knepley } 2022baeeceeSMatthew G. Knepley ((void (**)()) & wf->funcs->array[chunk.start])[ind] = func; 2036528b96dSMatthew G. Knepley PetscFunctionReturn(0); 2046528b96dSMatthew G. Knepley } 2056528b96dSMatthew G. Knepley 2069371c9d4SSatish Balay PetscErrorCode PetscWeakFormClearIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind) { 20706ad1575SMatthew G. Knepley PetscFormKey key; 20806ad1575SMatthew G. Knepley PetscChunk chunk; 20906ad1575SMatthew G. Knepley 21006ad1575SMatthew G. Knepley PetscFunctionBegin; 2119371c9d4SSatish Balay key.label = label; 2129371c9d4SSatish Balay key.value = value; 2139371c9d4SSatish Balay key.field = f; 2149371c9d4SSatish Balay key.part = part; 2159566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGet(ht, key, &chunk)); 21606ad1575SMatthew G. Knepley if (chunk.size < 0) { 21706ad1575SMatthew G. Knepley PetscFunctionReturn(0); 21806ad1575SMatthew G. Knepley } else if (!ind && chunk.size == 1) { 2199566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDel(ht, key)); 22006ad1575SMatthew G. Knepley PetscFunctionReturn(0); 22106ad1575SMatthew G. Knepley } else if (chunk.size <= ind) { 22206ad1575SMatthew G. Knepley PetscFunctionReturn(0); 22306ad1575SMatthew G. Knepley } 22406ad1575SMatthew G. Knepley ((void (**)()) & wf->funcs->array[chunk.start])[ind] = NULL; 22506ad1575SMatthew G. Knepley PetscFunctionReturn(0); 22606ad1575SMatthew G. Knepley } 22706ad1575SMatthew G. Knepley 22845480ffeSMatthew G. Knepley /*@ 22945480ffeSMatthew G. Knepley PetscWeakFormCopy - Copy the pointwise functions to another PetscWeakForm 23045480ffeSMatthew G. Knepley 23145480ffeSMatthew G. Knepley Not Collective 23245480ffeSMatthew G. Knepley 23345480ffeSMatthew G. Knepley Input Parameter: 23445480ffeSMatthew G. Knepley . wf - The original PetscWeakForm 23545480ffeSMatthew G. Knepley 23645480ffeSMatthew G. Knepley Output Parameter: 23745480ffeSMatthew G. Knepley . wfNew - The copy PetscWeakForm 23845480ffeSMatthew G. Knepley 23945480ffeSMatthew G. Knepley Level: intermediate 24045480ffeSMatthew G. Knepley 241db781477SPatrick Sanan .seealso: `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 24245480ffeSMatthew G. Knepley @*/ 2439371c9d4SSatish Balay PetscErrorCode PetscWeakFormCopy(PetscWeakForm wf, PetscWeakForm wfNew) { 24406ad1575SMatthew G. Knepley PetscInt f; 24545480ffeSMatthew G. Knepley 24645480ffeSMatthew G. Knepley PetscFunctionBegin; 24745480ffeSMatthew G. Knepley wfNew->Nf = wf->Nf; 2489566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferDestroy(&wfNew->funcs)); 2499566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferDuplicate(wf->funcs, &wfNew->funcs)); 25006ad1575SMatthew G. Knepley for (f = 0; f < PETSC_NUM_WF; ++f) { 2519566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDestroy(&wfNew->form[f])); 2529566063dSJacob Faibussowitsch PetscCall(PetscHMapFormDuplicate(wf->form[f], &wfNew->form[f])); 25306ad1575SMatthew G. Knepley } 25406ad1575SMatthew G. Knepley PetscFunctionReturn(0); 25506ad1575SMatthew G. Knepley } 25606ad1575SMatthew G. Knepley 25706ad1575SMatthew G. Knepley /*@ 25806ad1575SMatthew G. Knepley PetscWeakFormClear - Clear all functions from the PetscWeakForm 25906ad1575SMatthew G. Knepley 26006ad1575SMatthew G. Knepley Not Collective 26106ad1575SMatthew G. Knepley 26206ad1575SMatthew G. Knepley Input Parameter: 26306ad1575SMatthew G. Knepley . wf - The original PetscWeakForm 26406ad1575SMatthew G. Knepley 26506ad1575SMatthew G. Knepley Level: intermediate 26606ad1575SMatthew G. Knepley 267db781477SPatrick Sanan .seealso: `PetscWeakFormCopy()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 26806ad1575SMatthew G. Knepley @*/ 2699371c9d4SSatish Balay PetscErrorCode PetscWeakFormClear(PetscWeakForm wf) { 27006ad1575SMatthew G. Knepley PetscInt f; 27106ad1575SMatthew G. Knepley 27206ad1575SMatthew G. Knepley PetscFunctionBegin; 2739566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormClear(wf->form[f])); 27445480ffeSMatthew G. Knepley PetscFunctionReturn(0); 27545480ffeSMatthew G. Knepley } 27645480ffeSMatthew G. Knepley 2779371c9d4SSatish Balay static PetscErrorCode PetscWeakFormRewriteKeys_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label, PetscInt Nv, const PetscInt values[]) { 27806ad1575SMatthew G. Knepley PetscFormKey *keys; 27945480ffeSMatthew G. Knepley PetscInt n, i, v, off = 0; 28045480ffeSMatthew G. Knepley 28145480ffeSMatthew G. Knepley PetscFunctionBegin; 2829566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(hmap, &n)); 2839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &keys)); 2849566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetKeys(hmap, &off, keys)); 28545480ffeSMatthew G. Knepley for (i = 0; i < n; ++i) { 28645480ffeSMatthew G. Knepley if (keys[i].label == label) { 28745480ffeSMatthew G. Knepley PetscBool clear = PETSC_TRUE; 28845480ffeSMatthew G. Knepley void (**funcs)(); 28945480ffeSMatthew G. Knepley PetscInt Nf; 29045480ffeSMatthew G. Knepley 2919566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 29245480ffeSMatthew G. Knepley for (v = 0; v < Nv; ++v) { 2939566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, values[v], keys[i].field, keys[i].part, Nf, funcs)); 29445480ffeSMatthew G. Knepley if (values[v] == keys[i].value) clear = PETSC_FALSE; 29545480ffeSMatthew G. Knepley } 2969566063dSJacob Faibussowitsch if (clear) PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL)); 29745480ffeSMatthew G. Knepley } 29845480ffeSMatthew G. Knepley } 2999566063dSJacob Faibussowitsch PetscCall(PetscFree(keys)); 30045480ffeSMatthew G. Knepley PetscFunctionReturn(0); 30145480ffeSMatthew G. Knepley } 30245480ffeSMatthew G. Knepley 30345480ffeSMatthew G. Knepley /*@C 30445480ffeSMatthew G. Knepley PetscWeakFormRewriteKeys - Change any key on the given label to use the new set of label values 30545480ffeSMatthew G. Knepley 30645480ffeSMatthew G. Knepley Not Collective 30745480ffeSMatthew G. Knepley 30845480ffeSMatthew G. Knepley Input Parameters: 30945480ffeSMatthew G. Knepley + wf - The original PetscWeakForm 31045480ffeSMatthew G. Knepley . label - The label to change keys for 31145480ffeSMatthew G. Knepley . Nv - The number of new label values 31245480ffeSMatthew G. Knepley - values - The set of new values to relabel keys with 31345480ffeSMatthew G. Knepley 31445480ffeSMatthew G. Knepley Note: This is used internally when boundary label values are specified from the command line. 31545480ffeSMatthew G. Knepley 31645480ffeSMatthew G. Knepley Level: intermediate 31745480ffeSMatthew G. Knepley 318db781477SPatrick Sanan .seealso: `PetscWeakFormReplaceLabel()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 31945480ffeSMatthew G. Knepley @*/ 3209371c9d4SSatish Balay PetscErrorCode PetscWeakFormRewriteKeys(PetscWeakForm wf, DMLabel label, PetscInt Nv, const PetscInt values[]) { 32106ad1575SMatthew G. Knepley PetscInt f; 32245480ffeSMatthew G. Knepley 32345480ffeSMatthew G. Knepley PetscFunctionBegin; 3249566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormRewriteKeys_Internal(wf, wf->form[f], label, Nv, values)); 32545480ffeSMatthew G. Knepley PetscFunctionReturn(0); 32645480ffeSMatthew G. Knepley } 32745480ffeSMatthew G. Knepley 3289371c9d4SSatish Balay static PetscErrorCode PetscWeakFormReplaceLabel_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label) { 329d700741fSMatthew G. Knepley PetscFormKey *keys; 330d700741fSMatthew G. Knepley PetscInt n, i, off = 0, maxFuncs = 0; 331d700741fSMatthew G. Knepley void (**tmpf)(); 332d700741fSMatthew G. Knepley const char *name = NULL; 333d700741fSMatthew G. Knepley 334d700741fSMatthew G. Knepley PetscFunctionBegin; 3359566063dSJacob Faibussowitsch if (label) PetscCall(PetscObjectGetName((PetscObject)label, &name)); 3369566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(hmap, &n)); 3379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &keys)); 3389566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetKeys(hmap, &off, keys)); 339d700741fSMatthew G. Knepley for (i = 0; i < n; ++i) { 340d700741fSMatthew G. Knepley PetscBool match = PETSC_FALSE; 341d700741fSMatthew G. Knepley const char *lname = NULL; 342d700741fSMatthew G. Knepley 343d700741fSMatthew G. Knepley if (label == keys[i].label) continue; 3449566063dSJacob Faibussowitsch if (keys[i].label) PetscCall(PetscObjectGetName((PetscObject)keys[i].label, &lname)); 3459566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(name, lname, &match)); 346d700741fSMatthew G. Knepley if ((!name && !lname) || match) { 347d700741fSMatthew G. Knepley void (**funcs)(); 348d700741fSMatthew G. Knepley PetscInt Nf; 349d700741fSMatthew G. Knepley 3509566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs)); 351d700741fSMatthew G. Knepley maxFuncs = PetscMax(maxFuncs, Nf); 352d700741fSMatthew G. Knepley } 353d700741fSMatthew G. Knepley } 354d700741fSMatthew G. Knepley /* Need temp space because chunk buffer can be reallocated in SetFunction() call */ 3559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxFuncs, &tmpf)); 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, j; 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 for (j = 0; j < Nf; ++j) tmpf[j] = funcs[j]; 3699566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, label, keys[i].value, keys[i].field, keys[i].part, Nf, tmpf)); 3709566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL)); 371d700741fSMatthew G. Knepley } 372d700741fSMatthew G. Knepley } 3739566063dSJacob Faibussowitsch PetscCall(PetscFree(tmpf)); 3749566063dSJacob Faibussowitsch PetscCall(PetscFree(keys)); 375d700741fSMatthew G. Knepley PetscFunctionReturn(0); 376d700741fSMatthew G. Knepley } 377d700741fSMatthew G. Knepley 378d700741fSMatthew G. Knepley /*@C 379d700741fSMatthew G. Knepley PetscWeakFormReplaceLabel - Change any key on a label of the same name to use the new label 380d700741fSMatthew G. Knepley 381d700741fSMatthew G. Knepley Not Collective 382d700741fSMatthew G. Knepley 383d700741fSMatthew G. Knepley Input Parameters: 384d700741fSMatthew G. Knepley + wf - The original PetscWeakForm 385d700741fSMatthew G. Knepley - label - The label to change keys for 386d700741fSMatthew G. Knepley 387d700741fSMatthew G. Knepley Note: This is used internally when meshes are modified 388d700741fSMatthew G. Knepley 389d700741fSMatthew G. Knepley Level: intermediate 390d700741fSMatthew G. Knepley 391db781477SPatrick Sanan .seealso: `PetscWeakFormRewriteKeys()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()` 392d700741fSMatthew G. Knepley @*/ 3939371c9d4SSatish Balay PetscErrorCode PetscWeakFormReplaceLabel(PetscWeakForm wf, DMLabel label) { 394d700741fSMatthew G. Knepley PetscInt f; 395d700741fSMatthew G. Knepley 396d700741fSMatthew G. Knepley PetscFunctionBegin; 3979566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormReplaceLabel_Internal(wf, wf->form[f], label)); 398d700741fSMatthew G. Knepley PetscFunctionReturn(0); 399d700741fSMatthew G. Knepley } 400d700741fSMatthew G. Knepley 4019371c9d4SSatish Balay PetscErrorCode PetscWeakFormClearIndex(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscWeakFormKind kind, PetscInt ind) { 40206ad1575SMatthew G. Knepley PetscFunctionBegin; 4039566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClearIndexFunction_Private(wf, wf->form[kind], label, val, f, part, ind)); 40406ad1575SMatthew G. Knepley PetscFunctionReturn(0); 40506ad1575SMatthew G. Knepley } 40606ad1575SMatthew G. Knepley 4079371c9d4SSatish Balay 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[])) { 4086528b96dSMatthew G. Knepley PetscFunctionBegin; 4099566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (***)(void))obj)); 4106528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4116528b96dSMatthew G. Knepley } 4126528b96dSMatthew G. Knepley 4139371c9d4SSatish Balay 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[])) { 4146528b96dSMatthew G. Knepley PetscFunctionBegin; 4159566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (**)(void))obj)); 4166528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4176528b96dSMatthew G. Knepley } 4186528b96dSMatthew G. Knepley 4199371c9d4SSatish Balay 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[])) { 4206528b96dSMatthew G. Knepley PetscFunctionBegin; 4219566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, (void (*)(void))obj)); 4226528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4236528b96dSMatthew G. Knepley } 4246528b96dSMatthew G. Knepley 4259371c9d4SSatish Balay 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[])) { 4266528b96dSMatthew G. Knepley PetscFunctionBegin; 4279566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (void (**)(void))obj)); 4286528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4296528b96dSMatthew G. Knepley } 4306528b96dSMatthew G. Knepley 4319371c9d4SSatish Balay 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[])) { 4326528b96dSMatthew G. Knepley PetscFunctionBegin; 4339566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (void (*)(void))obj)); 4346528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4356528b96dSMatthew G. Knepley } 4366528b96dSMatthew G. Knepley 4379371c9d4SSatish Balay 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[])) { 4386528b96dSMatthew G. Knepley PetscFunctionBegin; 4399566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (***)(void))f0)); 4409566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (***)(void))f1)); 4416528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4426528b96dSMatthew G. Knepley } 4436528b96dSMatthew G. Knepley 4449371c9d4SSatish Balay 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[])) { 4456528b96dSMatthew G. Knepley PetscFunctionBegin; 4469566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, (void (*)(void))f0)); 4479566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, (void (*)(void))f1)); 4486528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4496528b96dSMatthew G. Knepley } 4506528b96dSMatthew G. Knepley 4519371c9d4SSatish Balay 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[])) { 4526528b96dSMatthew G. Knepley PetscFunctionBegin; 4539566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (**)(void))f0)); 4549566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (**)(void))f1)); 4556528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4566528b96dSMatthew G. Knepley } 4576528b96dSMatthew G. Knepley 4589371c9d4SSatish Balay 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[])) { 4596528b96dSMatthew G. Knepley PetscFunctionBegin; 4609566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, i0, (void (*)(void))f0)); 4619566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, i1, (void (*)(void))f1)); 4626528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4636528b96dSMatthew G. Knepley } 4646528b96dSMatthew G. Knepley 4659371c9d4SSatish Balay 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[])) { 4666528b96dSMatthew G. Knepley PetscFunctionBegin; 4679566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (***)(void))f0)); 4689566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (***)(void))f1)); 4696528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4706528b96dSMatthew G. Knepley } 4716528b96dSMatthew G. Knepley 4729371c9d4SSatish Balay 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[])) { 4736528b96dSMatthew G. Knepley PetscFunctionBegin; 4749566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, (void (*)(void))f0)); 4759566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, (void (*)(void))f1)); 4766528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4776528b96dSMatthew G. Knepley } 4786528b96dSMatthew G. Knepley 4799371c9d4SSatish Balay 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[])) { 4806528b96dSMatthew G. Knepley PetscFunctionBegin; 4819566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (**)(void))f0)); 4829566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (**)(void))f1)); 4836528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4846528b96dSMatthew G. Knepley } 4856528b96dSMatthew G. Knepley 4869371c9d4SSatish Balay 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[])) { 4876528b96dSMatthew G. Knepley PetscFunctionBegin; 4889566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, i0, (void (*)(void))f0)); 4899566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, i1, (void (*)(void))f1)); 4906528b96dSMatthew G. Knepley PetscFunctionReturn(0); 4916528b96dSMatthew G. Knepley } 4926528b96dSMatthew G. Knepley 4939371c9d4SSatish Balay PetscErrorCode PetscWeakFormHasJacobian(PetscWeakForm wf, PetscBool *hasJac) { 4946528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 4956528b96dSMatthew G. Knepley 4966528b96dSMatthew G. Knepley PetscFunctionBegin; 4976528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 4986528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJac, 2); 4999566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G0], &n0)); 5009566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G1], &n1)); 5019566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G2], &n2)); 5029566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G3], &n3)); 5036528b96dSMatthew G. Knepley *hasJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 5046528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5056528b96dSMatthew G. Knepley } 5066528b96dSMatthew G. Knepley 5079371c9d4SSatish Balay 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[])) { 5086528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5096528b96dSMatthew G. Knepley 5106528b96dSMatthew G. Knepley PetscFunctionBegin; 5119566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (***)(void))g0)); 5129566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (***)(void))g1)); 5139566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (***)(void))g2)); 5149566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (***)(void))g3)); 5156528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5166528b96dSMatthew G. Knepley } 5176528b96dSMatthew G. Knepley 5189371c9d4SSatish Balay 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[])) { 5196528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5206528b96dSMatthew G. Knepley 5216528b96dSMatthew G. Knepley PetscFunctionBegin; 5229566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, (void (*)(void))g0)); 5239566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, (void (*)(void))g1)); 5249566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, (void (*)(void))g2)); 5259566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, (void (*)(void))g3)); 5266528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5276528b96dSMatthew G. Knepley } 5286528b96dSMatthew G. Knepley 5299371c9d4SSatish Balay 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[])) { 5306528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5316528b96dSMatthew G. Knepley 5326528b96dSMatthew G. Knepley PetscFunctionBegin; 5339566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (**)(void))g0)); 5349566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (**)(void))g1)); 5359566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (**)(void))g2)); 5369566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (**)(void))g3)); 5376528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5386528b96dSMatthew G. Knepley } 5396528b96dSMatthew G. Knepley 5409371c9d4SSatish Balay 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[])) { 5416528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5426528b96dSMatthew G. Knepley 5436528b96dSMatthew G. Knepley PetscFunctionBegin; 5449566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, i0, (void (*)(void))g0)); 5459566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, i1, (void (*)(void))g1)); 5469566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, i2, (void (*)(void))g2)); 5479566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, i3, (void (*)(void))g3)); 5486528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5496528b96dSMatthew G. Knepley } 5506528b96dSMatthew G. Knepley 5519371c9d4SSatish Balay PetscErrorCode PetscWeakFormHasJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre) { 5526528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 5536528b96dSMatthew G. Knepley 5546528b96dSMatthew G. Knepley PetscFunctionBegin; 5556528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 5566528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJacPre, 2); 5579566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP0], &n0)); 5589566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP1], &n1)); 5599566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP2], &n2)); 5609566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP3], &n3)); 5616528b96dSMatthew G. Knepley *hasJacPre = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 5626528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5636528b96dSMatthew G. Knepley } 5646528b96dSMatthew G. Knepley 5659371c9d4SSatish Balay 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[])) { 5666528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5676528b96dSMatthew G. Knepley 5686528b96dSMatthew G. Knepley PetscFunctionBegin; 5699566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (***)(void))g0)); 5709566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (***)(void))g1)); 5719566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (***)(void))g2)); 5729566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (***)(void))g3)); 5736528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5746528b96dSMatthew G. Knepley } 5756528b96dSMatthew G. Knepley 5769371c9d4SSatish Balay 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[])) { 5776528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5786528b96dSMatthew G. Knepley 5796528b96dSMatthew G. Knepley PetscFunctionBegin; 5809566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, (void (*)(void))g0)); 5819566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, (void (*)(void))g1)); 5829566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, (void (*)(void))g2)); 5839566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, (void (*)(void))g3)); 5846528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5856528b96dSMatthew G. Knepley } 5866528b96dSMatthew G. Knepley 5879371c9d4SSatish Balay 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[])) { 5886528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 5896528b96dSMatthew G. Knepley 5906528b96dSMatthew G. Knepley PetscFunctionBegin; 5919566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (**)(void))g0)); 5929566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (**)(void))g1)); 5939566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (**)(void))g2)); 5949566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (**)(void))g3)); 5956528b96dSMatthew G. Knepley PetscFunctionReturn(0); 5966528b96dSMatthew G. Knepley } 5976528b96dSMatthew G. Knepley 5989371c9d4SSatish Balay 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[])) { 5996528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6006528b96dSMatthew G. Knepley 6016528b96dSMatthew G. Knepley PetscFunctionBegin; 6029566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, i0, (void (*)(void))g0)); 6039566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, i1, (void (*)(void))g1)); 6049566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, i2, (void (*)(void))g2)); 6059566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, i3, (void (*)(void))g3)); 6066528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6076528b96dSMatthew G. Knepley } 6086528b96dSMatthew G. Knepley 6099371c9d4SSatish Balay PetscErrorCode PetscWeakFormHasBdJacobian(PetscWeakForm wf, PetscBool *hasJac) { 6106528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 6116528b96dSMatthew G. Knepley 6126528b96dSMatthew G. Knepley PetscFunctionBegin; 6136528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 6146528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJac, 2); 6159566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG0], &n0)); 6169566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG1], &n1)); 6179566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG2], &n2)); 6189566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG3], &n3)); 6196528b96dSMatthew G. Knepley *hasJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 6206528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6216528b96dSMatthew G. Knepley } 6226528b96dSMatthew G. Knepley 6239371c9d4SSatish Balay 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[])) { 6246528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6256528b96dSMatthew G. Knepley 6266528b96dSMatthew G. Knepley PetscFunctionBegin; 6279566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (***)(void))g0)); 6289566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (***)(void))g1)); 6299566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (***)(void))g2)); 6309566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (***)(void))g3)); 6316528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6326528b96dSMatthew G. Knepley } 6336528b96dSMatthew G. Knepley 6349371c9d4SSatish Balay 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[])) { 6356528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6366528b96dSMatthew G. Knepley 6376528b96dSMatthew G. Knepley PetscFunctionBegin; 6389566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, (void (*)(void))g0)); 6399566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, (void (*)(void))g1)); 6409566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, (void (*)(void))g2)); 6419566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, (void (*)(void))g3)); 6426528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6436528b96dSMatthew G. Knepley } 6446528b96dSMatthew G. Knepley 6459371c9d4SSatish Balay 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[])) { 6466528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6476528b96dSMatthew G. Knepley 6486528b96dSMatthew G. Knepley PetscFunctionBegin; 6499566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (**)(void))g0)); 6509566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (**)(void))g1)); 6519566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (**)(void))g2)); 6529566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (**)(void))g3)); 6536528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6546528b96dSMatthew G. Knepley } 6556528b96dSMatthew G. Knepley 6569371c9d4SSatish Balay 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[])) { 6576528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6586528b96dSMatthew G. Knepley 6596528b96dSMatthew G. Knepley PetscFunctionBegin; 6609566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, i0, (void (*)(void))g0)); 6619566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, i1, (void (*)(void))g1)); 6629566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, i2, (void (*)(void))g2)); 6639566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, i3, (void (*)(void))g3)); 6646528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6656528b96dSMatthew G. Knepley } 6666528b96dSMatthew G. Knepley 6679371c9d4SSatish Balay PetscErrorCode PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre) { 6686528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 6696528b96dSMatthew G. Knepley 6706528b96dSMatthew G. Knepley PetscFunctionBegin; 6716528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 6726528b96dSMatthew G. Knepley PetscValidBoolPointer(hasJacPre, 2); 6739566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP0], &n0)); 6749566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP1], &n1)); 6759566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP2], &n2)); 6769566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP3], &n3)); 6776528b96dSMatthew G. Knepley *hasJacPre = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 6786528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6796528b96dSMatthew G. Knepley } 6806528b96dSMatthew G. Knepley 6819371c9d4SSatish Balay 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[])) { 6826528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6836528b96dSMatthew G. Knepley 6846528b96dSMatthew G. Knepley PetscFunctionBegin; 6859566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (***)(void))g0)); 6869566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (***)(void))g1)); 6879566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (***)(void))g2)); 6889566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (***)(void))g3)); 6896528b96dSMatthew G. Knepley PetscFunctionReturn(0); 6906528b96dSMatthew G. Knepley } 6916528b96dSMatthew G. Knepley 6929371c9d4SSatish Balay 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[])) { 6936528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 6946528b96dSMatthew G. Knepley 6956528b96dSMatthew G. Knepley PetscFunctionBegin; 6969566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, (void (*)(void))g0)); 6979566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, (void (*)(void))g1)); 6989566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, (void (*)(void))g2)); 6999566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, (void (*)(void))g3)); 7006528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7016528b96dSMatthew G. Knepley } 7026528b96dSMatthew G. Knepley 7039371c9d4SSatish Balay 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[])) { 7046528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7056528b96dSMatthew G. Knepley 7066528b96dSMatthew G. Knepley PetscFunctionBegin; 7079566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (**)(void))g0)); 7089566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (**)(void))g1)); 7099566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (**)(void))g2)); 7109566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (**)(void))g3)); 7116528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7126528b96dSMatthew G. Knepley } 7136528b96dSMatthew G. Knepley 7149371c9d4SSatish Balay 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[])) { 7156528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7166528b96dSMatthew G. Knepley 7176528b96dSMatthew G. Knepley PetscFunctionBegin; 7189566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, i0, (void (*)(void))g0)); 7199566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, i1, (void (*)(void))g1)); 7209566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, i2, (void (*)(void))g2)); 7219566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, i3, (void (*)(void))g3)); 7226528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7236528b96dSMatthew G. Knepley } 7246528b96dSMatthew G. Knepley 7259371c9d4SSatish Balay PetscErrorCode PetscWeakFormHasDynamicJacobian(PetscWeakForm wf, PetscBool *hasDynJac) { 7266528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3; 7276528b96dSMatthew G. Knepley 7286528b96dSMatthew G. Knepley PetscFunctionBegin; 7296528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 7306528b96dSMatthew G. Knepley PetscValidBoolPointer(hasDynJac, 2); 7319566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT0], &n0)); 7329566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT1], &n1)); 7339566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT2], &n2)); 7349566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT3], &n3)); 7356528b96dSMatthew G. Knepley *hasDynJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE; 7366528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7376528b96dSMatthew G. Knepley } 7386528b96dSMatthew G. Knepley 7399371c9d4SSatish Balay 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[])) { 7406528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7416528b96dSMatthew G. Knepley 7426528b96dSMatthew G. Knepley PetscFunctionBegin; 7439566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (***)(void))g0)); 7449566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (***)(void))g1)); 7459566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (***)(void))g2)); 7469566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (***)(void))g3)); 7476528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7486528b96dSMatthew G. Knepley } 7496528b96dSMatthew G. Knepley 7509371c9d4SSatish Balay 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[])) { 7516528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7526528b96dSMatthew G. Knepley 7536528b96dSMatthew G. Knepley PetscFunctionBegin; 7549566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, (void (*)(void))g0)); 7559566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, (void (*)(void))g1)); 7569566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, (void (*)(void))g2)); 7579566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, (void (*)(void))g3)); 7586528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7596528b96dSMatthew G. Knepley } 7606528b96dSMatthew G. Knepley 7619371c9d4SSatish Balay 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[])) { 7626528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7636528b96dSMatthew G. Knepley 7646528b96dSMatthew G. Knepley PetscFunctionBegin; 7659566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (**)(void))g0)); 7669566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (**)(void))g1)); 7679566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (**)(void))g2)); 7689566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (**)(void))g3)); 7696528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7706528b96dSMatthew G. Knepley } 7716528b96dSMatthew G. Knepley 7729371c9d4SSatish Balay 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[])) { 7736528b96dSMatthew G. Knepley PetscInt find = f * wf->Nf + g; 7746528b96dSMatthew G. Knepley 7756528b96dSMatthew G. Knepley PetscFunctionBegin; 7769566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, i0, (void (*)(void))g0)); 7779566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, i1, (void (*)(void))g1)); 7789566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, i2, (void (*)(void))g2)); 7799566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, i3, (void (*)(void))g3)); 7806528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7816528b96dSMatthew G. Knepley } 7826528b96dSMatthew G. Knepley 7839371c9d4SSatish Balay 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 *)) { 7846528b96dSMatthew G. Knepley PetscFunctionBegin; 7859566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (***)(void))r)); 7866528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7876528b96dSMatthew G. Knepley } 7886528b96dSMatthew G. Knepley 7899371c9d4SSatish Balay 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 *)) { 7906528b96dSMatthew G. Knepley PetscFunctionBegin; 7919566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (**)(void))r)); 7926528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7936528b96dSMatthew G. Knepley } 7946528b96dSMatthew G. Knepley 7959371c9d4SSatish Balay 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 *)) { 7966528b96dSMatthew G. Knepley PetscFunctionBegin; 7979566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, i, (void (*)(void))r)); 7986528b96dSMatthew G. Knepley PetscFunctionReturn(0); 7996528b96dSMatthew G. Knepley } 8006528b96dSMatthew G. Knepley 8016528b96dSMatthew G. Knepley /*@ 8026528b96dSMatthew G. Knepley PetscWeakFormGetNumFields - Returns the number of fields 8036528b96dSMatthew G. Knepley 8046528b96dSMatthew G. Knepley Not collective 8056528b96dSMatthew G. Knepley 8066528b96dSMatthew G. Knepley Input Parameter: 8076528b96dSMatthew G. Knepley . wf - The PetscWeakForm object 8086528b96dSMatthew G. Knepley 8096528b96dSMatthew G. Knepley Output Parameter: 810a5b23f4aSJose E. Roman . Nf - The number of fields 8116528b96dSMatthew G. Knepley 8126528b96dSMatthew G. Knepley Level: beginner 8136528b96dSMatthew G. Knepley 814db781477SPatrick Sanan .seealso: `PetscWeakFormSetNumFields()`, `PetscWeakFormCreate()` 8156528b96dSMatthew G. Knepley @*/ 8169371c9d4SSatish Balay PetscErrorCode PetscWeakFormGetNumFields(PetscWeakForm wf, PetscInt *Nf) { 8176528b96dSMatthew G. Knepley PetscFunctionBegin; 8186528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 819dadcf809SJacob Faibussowitsch PetscValidIntPointer(Nf, 2); 8206528b96dSMatthew G. Knepley *Nf = wf->Nf; 8216528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8226528b96dSMatthew G. Knepley } 8236528b96dSMatthew G. Knepley 8246528b96dSMatthew G. Knepley /*@ 8256528b96dSMatthew G. Knepley PetscWeakFormSetNumFields - Sets the number of fields 8266528b96dSMatthew G. Knepley 8276528b96dSMatthew G. Knepley Not collective 8286528b96dSMatthew G. Knepley 8296528b96dSMatthew G. Knepley Input Parameters: 8306528b96dSMatthew G. Knepley + wf - The PetscWeakForm object 8316528b96dSMatthew G. Knepley - Nf - The number of fields 8326528b96dSMatthew G. Knepley 8336528b96dSMatthew G. Knepley Level: beginner 8346528b96dSMatthew G. Knepley 835db781477SPatrick Sanan .seealso: `PetscWeakFormGetNumFields()`, `PetscWeakFormCreate()` 8366528b96dSMatthew G. Knepley @*/ 8379371c9d4SSatish Balay PetscErrorCode PetscWeakFormSetNumFields(PetscWeakForm wf, PetscInt Nf) { 8386528b96dSMatthew G. Knepley PetscFunctionBegin; 8396528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 8406528b96dSMatthew G. Knepley wf->Nf = Nf; 8416528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8426528b96dSMatthew G. Knepley } 8436528b96dSMatthew G. Knepley 8446528b96dSMatthew G. Knepley /*@ 8456528b96dSMatthew G. Knepley PetscWeakFormDestroy - Destroys a PetscWeakForm object 8466528b96dSMatthew G. Knepley 8476528b96dSMatthew G. Knepley Collective on wf 8486528b96dSMatthew G. Knepley 8496528b96dSMatthew G. Knepley Input Parameter: 8506528b96dSMatthew G. Knepley . wf - the PetscWeakForm object to destroy 8516528b96dSMatthew G. Knepley 8526528b96dSMatthew G. Knepley Level: developer 8536528b96dSMatthew G. Knepley 854db781477SPatrick Sanan .seealso `PetscWeakFormCreate()`, `PetscWeakFormView()` 8556528b96dSMatthew G. Knepley @*/ 8569371c9d4SSatish Balay PetscErrorCode PetscWeakFormDestroy(PetscWeakForm *wf) { 85706ad1575SMatthew G. Knepley PetscInt f; 8586528b96dSMatthew G. Knepley 8596528b96dSMatthew G. Knepley PetscFunctionBegin; 8606528b96dSMatthew G. Knepley if (!*wf) PetscFunctionReturn(0); 8616528b96dSMatthew G. Knepley PetscValidHeaderSpecific((*wf), PETSCWEAKFORM_CLASSID, 1); 8626528b96dSMatthew G. Knepley 8639371c9d4SSatish Balay if (--((PetscObject)(*wf))->refct > 0) { 8649371c9d4SSatish Balay *wf = NULL; 8659371c9d4SSatish Balay PetscFunctionReturn(0); 8669371c9d4SSatish Balay } 8676528b96dSMatthew G. Knepley ((PetscObject)(*wf))->refct = 0; 8689566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferDestroy(&(*wf)->funcs)); 8699566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormDestroy(&(*wf)->form[f])); 8709566063dSJacob Faibussowitsch PetscCall(PetscFree((*wf)->form)); 8719566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(wf)); 8726528b96dSMatthew G. Knepley PetscFunctionReturn(0); 8736528b96dSMatthew G. Knepley } 8746528b96dSMatthew G. Knepley 8759371c9d4SSatish Balay static PetscErrorCode PetscWeakFormViewTable_Ascii(PetscWeakForm wf, PetscViewer viewer, PetscBool splitField, const char tableName[], PetscHMapForm map) { 87645480ffeSMatthew G. Knepley PetscInt Nf = wf->Nf, Nk, k; 8776528b96dSMatthew G. Knepley 8786528b96dSMatthew G. Knepley PetscFunctionBegin; 8799566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(map, &Nk)); 8806528b96dSMatthew G. Knepley if (Nk) { 88106ad1575SMatthew G. Knepley PetscFormKey *keys; 882165f9cc3SJed Brown void (**funcs)(void) = NULL; 8835fedec97SMatthew G. Knepley const char **names; 8845fedec97SMatthew G. Knepley PetscInt *values, *idx1, *idx2, *idx; 8855fedec97SMatthew G. Knepley PetscBool showPart = PETSC_FALSE, showPointer = PETSC_FALSE; 8865fedec97SMatthew G. Knepley PetscInt off = 0; 8876528b96dSMatthew G. Knepley 8889566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(Nk, &keys, Nk, &names, Nk, &values, Nk, &idx1, Nk, &idx2, Nk, &idx)); 8899566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetKeys(map, &off, keys)); 8905fedec97SMatthew G. Knepley /* Sort keys by label name and value */ 8915fedec97SMatthew G. Knepley { 8925fedec97SMatthew G. Knepley /* First sort values */ 8939371c9d4SSatish Balay for (k = 0; k < Nk; ++k) { 8949371c9d4SSatish Balay values[k] = keys[k].value; 8959371c9d4SSatish Balay idx1[k] = k; 8969371c9d4SSatish Balay } 8979566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(Nk, values, idx1)); 8985fedec97SMatthew G. Knepley /* If the string sort is stable, it will be sorted correctly overall */ 8995fedec97SMatthew G. Knepley for (k = 0; k < Nk; ++k) { 9009566063dSJacob Faibussowitsch if (keys[idx1[k]].label) PetscCall(PetscObjectGetName((PetscObject)keys[idx1[k]].label, &names[k])); 9015fedec97SMatthew G. Knepley else { names[k] = ""; } 9025fedec97SMatthew G. Knepley idx2[k] = k; 9035fedec97SMatthew G. Knepley } 9049566063dSJacob Faibussowitsch PetscCall(PetscSortStrWithPermutation(Nk, names, idx2)); 9055fedec97SMatthew G. Knepley for (k = 0; k < Nk; ++k) { 9069566063dSJacob Faibussowitsch if (keys[k].label) PetscCall(PetscObjectGetName((PetscObject)keys[k].label, &names[k])); 9075fedec97SMatthew G. Knepley else { names[k] = ""; } 9085fedec97SMatthew G. Knepley idx[k] = idx1[idx2[k]]; 9095fedec97SMatthew G. Knepley } 9105fedec97SMatthew G. Knepley } 9119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", tableName)); 9129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 9136528b96dSMatthew G. Knepley for (k = 0; k < Nk; ++k) { 9140944a4d6SMatthew G. Knepley if (keys[k].part != 0) showPart = PETSC_TRUE; 9150944a4d6SMatthew G. Knepley } 9160944a4d6SMatthew G. Knepley for (k = 0; k < Nk; ++k) { 9175fedec97SMatthew G. Knepley const PetscInt i = idx[k]; 9185fedec97SMatthew G. Knepley PetscInt n, f; 9195fedec97SMatthew G. Knepley 9205fedec97SMatthew G. Knepley if (keys[i].label) { 92163a3b9bcSJacob Faibussowitsch if (showPointer) PetscCall(PetscViewerASCIIPrintf(viewer, "(%s:%p, %" PetscInt_FMT ") ", names[i], keys[i].label, keys[i].value)); 92263a3b9bcSJacob Faibussowitsch else PetscCall(PetscViewerASCIIPrintf(viewer, "(%s, %" PetscInt_FMT ") ", names[i], keys[i].value)); 92363a3b9bcSJacob Faibussowitsch } 9249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 92563a3b9bcSJacob Faibussowitsch if (splitField) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ", %" PetscInt_FMT ") ", keys[i].field / Nf, keys[i].field % Nf)); 92663a3b9bcSJacob Faibussowitsch else PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") ", keys[i].field)); 92763a3b9bcSJacob Faibussowitsch if (showPart) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") ", keys[i].part)); 9289566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetFunction_Private(wf, map, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &n, &funcs)); 9295fedec97SMatthew G. Knepley for (f = 0; f < n; ++f) { 930258ec3d2SMatthew G. Knepley char *fname; 9315fedec97SMatthew G. Knepley size_t len, l; 932258ec3d2SMatthew G. Knepley 9339566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 9349566063dSJacob Faibussowitsch PetscCall(PetscDLAddr(funcs[f], &fname)); 9355fedec97SMatthew G. Knepley if (fname) { 9365fedec97SMatthew G. Knepley /* Eliminate argument types */ 9379566063dSJacob Faibussowitsch PetscCall(PetscStrlen(fname, &len)); 9389371c9d4SSatish Balay for (l = 0; l < len; ++l) 9399371c9d4SSatish Balay if (fname[l] == '(') { 9409371c9d4SSatish Balay fname[l] = '\0'; 9419371c9d4SSatish Balay break; 9429371c9d4SSatish Balay } 9439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%s", fname)); 9445fedec97SMatthew G. Knepley } else if (showPointer) { 9459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%p", funcs[f])); 9465fedec97SMatthew G. Knepley } 9479566063dSJacob Faibussowitsch PetscCall(PetscFree(fname)); 9486528b96dSMatthew G. Knepley } 9499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 9509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 9516528b96dSMatthew G. Knepley } 9529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 9539566063dSJacob Faibussowitsch PetscCall(PetscFree6(keys, names, values, idx1, idx2, idx)); 9546528b96dSMatthew G. Knepley } 9556528b96dSMatthew G. Knepley PetscFunctionReturn(0); 9566528b96dSMatthew G. Knepley } 9576528b96dSMatthew G. Knepley 9589371c9d4SSatish Balay static PetscErrorCode PetscWeakFormView_Ascii(PetscWeakForm wf, PetscViewer viewer) { 9596528b96dSMatthew G. Knepley PetscViewerFormat format; 96006ad1575SMatthew G. Knepley PetscInt f; 9616528b96dSMatthew G. Knepley 9626528b96dSMatthew G. Knepley PetscFunctionBegin; 9639566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 96463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Weak Form System with %" PetscInt_FMT " fields\n", wf->Nf)); 9659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 966*48a46eb9SPierre Jolivet for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormViewTable_Ascii(wf, viewer, PETSC_TRUE, PetscWeakFormKinds[f], wf->form[f])); 9679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 9686528b96dSMatthew G. Knepley PetscFunctionReturn(0); 9696528b96dSMatthew G. Knepley } 9706528b96dSMatthew G. Knepley 9716528b96dSMatthew G. Knepley /*@C 9726528b96dSMatthew G. Knepley PetscWeakFormView - Views a PetscWeakForm 9736528b96dSMatthew G. Knepley 9746528b96dSMatthew G. Knepley Collective on wf 9756528b96dSMatthew G. Knepley 976d8d19677SJose E. Roman Input Parameters: 9776528b96dSMatthew G. Knepley + wf - the PetscWeakForm object to view 9786528b96dSMatthew G. Knepley - v - the viewer 9796528b96dSMatthew G. Knepley 9806528b96dSMatthew G. Knepley Level: developer 9816528b96dSMatthew G. Knepley 982db781477SPatrick Sanan .seealso `PetscWeakFormDestroy()`, `PetscWeakFormCreate()` 9836528b96dSMatthew G. Knepley @*/ 9849371c9d4SSatish Balay PetscErrorCode PetscWeakFormView(PetscWeakForm wf, PetscViewer v) { 9856528b96dSMatthew G. Knepley PetscBool iascii; 9866528b96dSMatthew G. Knepley 9876528b96dSMatthew G. Knepley PetscFunctionBegin; 9886528b96dSMatthew G. Knepley PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1); 9899566063dSJacob Faibussowitsch if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)wf), &v)); 9906528b96dSMatthew G. Knepley else { PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2); } 9919566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii)); 9929566063dSJacob Faibussowitsch if (iascii) PetscCall(PetscWeakFormView_Ascii(wf, v)); 993dbbe0bcdSBarry Smith PetscTryTypeMethod(wf, view, v); 9946528b96dSMatthew G. Knepley PetscFunctionReturn(0); 9956528b96dSMatthew G. Knepley } 9966528b96dSMatthew G. Knepley 9976528b96dSMatthew G. Knepley /*@ 9986528b96dSMatthew G. Knepley PetscWeakFormCreate - Creates an empty PetscWeakForm object. 9996528b96dSMatthew G. Knepley 10006528b96dSMatthew G. Knepley Collective 10016528b96dSMatthew G. Knepley 10026528b96dSMatthew G. Knepley Input Parameter: 10036528b96dSMatthew G. Knepley . comm - The communicator for the PetscWeakForm object 10046528b96dSMatthew G. Knepley 10056528b96dSMatthew G. Knepley Output Parameter: 10066528b96dSMatthew G. Knepley . wf - The PetscWeakForm object 10076528b96dSMatthew G. Knepley 10086528b96dSMatthew G. Knepley Level: beginner 10096528b96dSMatthew G. Knepley 1010db781477SPatrick Sanan .seealso: `PetscDS`, `PetscWeakFormDestroy()` 10116528b96dSMatthew G. Knepley @*/ 10129371c9d4SSatish Balay PetscErrorCode PetscWeakFormCreate(MPI_Comm comm, PetscWeakForm *wf) { 10136528b96dSMatthew G. Knepley PetscWeakForm p; 101406ad1575SMatthew G. Knepley PetscInt f; 10156528b96dSMatthew G. Knepley 10166528b96dSMatthew G. Knepley PetscFunctionBegin; 10176528b96dSMatthew G. Knepley PetscValidPointer(wf, 2); 10186528b96dSMatthew G. Knepley *wf = NULL; 10199566063dSJacob Faibussowitsch PetscCall(PetscDSInitializePackage()); 10206528b96dSMatthew G. Knepley 10219566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(p, PETSCWEAKFORM_CLASSID, "PetscWeakForm", "Weak Form System", "PetscWeakForm", comm, PetscWeakFormDestroy, PetscWeakFormView)); 10226528b96dSMatthew G. Knepley 10236528b96dSMatthew G. Knepley p->Nf = 0; 10249566063dSJacob Faibussowitsch PetscCall(PetscChunkBufferCreate(sizeof(&PetscWeakFormCreate), 2, &p->funcs)); 10259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PETSC_NUM_WF, &p->form)); 10269566063dSJacob Faibussowitsch for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormCreate(&p->form[f])); 10276528b96dSMatthew G. Knepley *wf = p; 10286528b96dSMatthew G. Knepley PetscFunctionReturn(0); 10296528b96dSMatthew G. Knepley } 1030