106ad1575SMatthew G. Knepley static char help[] = "Tests for PetscWeakForm.\n\n";
206ad1575SMatthew G. Knepley
306ad1575SMatthew G. Knepley #include <petscds.h>
406ad1575SMatthew G. Knepley
f0(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])5d71ae5a4SJacob Faibussowitsch static void f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
6d71ae5a4SJacob Faibussowitsch {
706ad1575SMatthew G. Knepley f0[0] = 0.0;
806ad1575SMatthew G. Knepley }
906ad1575SMatthew G. Knepley
f1(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])10d71ae5a4SJacob Faibussowitsch static void f1(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
11d71ae5a4SJacob Faibussowitsch {
1206ad1575SMatthew G. Knepley f0[0] = 0.0;
1306ad1575SMatthew G. Knepley }
1406ad1575SMatthew G. Knepley
f2(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])15d71ae5a4SJacob Faibussowitsch static void f2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
16d71ae5a4SJacob Faibussowitsch {
1706ad1575SMatthew G. Knepley f0[0] = 0.0;
1806ad1575SMatthew G. Knepley }
1906ad1575SMatthew G. Knepley
f3(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])20d71ae5a4SJacob Faibussowitsch static void f3(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
21d71ae5a4SJacob Faibussowitsch {
2206ad1575SMatthew G. Knepley f0[0] = 0.0;
2306ad1575SMatthew G. Knepley }
2406ad1575SMatthew G. Knepley
CheckResidual(PetscWeakForm wf,PetscFormKey key,PetscInt in0,PetscPointFn * if0[],PetscInt in1,PetscPointFn * if1[])252192575eSBarry Smith static PetscErrorCode CheckResidual(PetscWeakForm wf, PetscFormKey key, PetscInt in0, PetscPointFn *if0[], PetscInt in1, PetscPointFn *if1[])
26d71ae5a4SJacob Faibussowitsch {
272192575eSBarry Smith PetscPointFn **f0, **f1;
2806ad1575SMatthew G. Knepley PetscInt n0, n1, i;
2906ad1575SMatthew G. Knepley
3006ad1575SMatthew G. Knepley PetscFunctionBegin;
319566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0, &n1, &f1));
3263a3b9bcSJacob Faibussowitsch PetscCheck(n0 == in0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Found %" PetscInt_FMT " f0 functions != %" PetscInt_FMT " functions input", n0, in0);
33ad540459SPierre Jolivet for (i = 0; i < n0; ++i) PetscCheck(f0[i] == if0[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "f0[%" PetscInt_FMT "] != input f0[%" PetscInt_FMT "]", i, i);
3463a3b9bcSJacob Faibussowitsch PetscCheck(n1 == in1, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Found %" PetscInt_FMT " f1 functions != %" PetscInt_FMT " functions input", n0, in0);
35ad540459SPierre Jolivet for (i = 0; i < n1; ++i) PetscCheck(f1[i] == if1[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "f1[%" PetscInt_FMT "] != input f1[%" PetscInt_FMT "]", i, i);
363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3706ad1575SMatthew G. Knepley }
3806ad1575SMatthew G. Knepley
TestSetIndex(PetscWeakForm wf)39d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestSetIndex(PetscWeakForm wf)
40d71ae5a4SJacob Faibussowitsch {
412192575eSBarry Smith PetscPointFn *f[4] = {f0, f1, f2, f3};
4206ad1575SMatthew G. Knepley DMLabel label;
4306ad1575SMatthew G. Knepley const PetscInt value = 3, field = 1, part = 2;
4406ad1575SMatthew G. Knepley PetscFormKey key;
4506ad1575SMatthew G. Knepley PetscInt i, j, k, l;
4606ad1575SMatthew G. Knepley
4706ad1575SMatthew G. Knepley PetscFunctionBegin;
489566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
499371c9d4SSatish Balay key.label = label;
509371c9d4SSatish Balay key.value = value;
519371c9d4SSatish Balay key.field = field;
529371c9d4SSatish Balay key.part = part;
5306ad1575SMatthew G. Knepley /* Check f0 */
5406ad1575SMatthew G. Knepley for (i = 0; i < 4; ++i) {
5506ad1575SMatthew G. Knepley for (j = 0; j < 4; ++j) {
5606ad1575SMatthew G. Knepley if (j == i) continue;
5706ad1575SMatthew G. Knepley for (k = 0; k < 4; ++k) {
5806ad1575SMatthew G. Knepley if ((k == i) || (k == j)) continue;
5906ad1575SMatthew G. Knepley for (l = 0; l < 4; ++l) {
6006ad1575SMatthew G. Knepley if ((l == i) || (l == j) || (l == k)) continue;
619566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], 0, NULL));
629566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], 0, NULL));
639566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], 0, NULL));
649566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], 0, NULL));
659566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
669566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
6706ad1575SMatthew G. Knepley }
6806ad1575SMatthew G. Knepley }
6906ad1575SMatthew G. Knepley }
7006ad1575SMatthew G. Knepley }
7106ad1575SMatthew G. Knepley /* Check f1 */
7206ad1575SMatthew G. Knepley for (i = 0; i < 4; ++i) {
7306ad1575SMatthew G. Knepley for (j = 0; j < 4; ++j) {
7406ad1575SMatthew G. Knepley if (j == i) continue;
7506ad1575SMatthew G. Knepley for (k = 0; k < 4; ++k) {
7606ad1575SMatthew G. Knepley if ((k == i) || (k == j)) continue;
7706ad1575SMatthew G. Knepley for (l = 0; l < 4; ++l) {
7806ad1575SMatthew G. Knepley if ((l == i) || (l == j) || (l == k)) continue;
799566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, i, f[i]));
809566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, j, f[j]));
819566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, k, f[k]));
829566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, l, f[l]));
839566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
849566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
8506ad1575SMatthew G. Knepley }
8606ad1575SMatthew G. Knepley }
8706ad1575SMatthew G. Knepley }
8806ad1575SMatthew G. Knepley }
8906ad1575SMatthew G. Knepley /* Check f0 and f1 */
9006ad1575SMatthew G. Knepley for (i = 0; i < 4; ++i) {
9106ad1575SMatthew G. Knepley for (j = 0; j < 4; ++j) {
9206ad1575SMatthew G. Knepley if (j == i) continue;
9306ad1575SMatthew G. Knepley for (k = 0; k < 4; ++k) {
9406ad1575SMatthew G. Knepley if ((k == i) || (k == j)) continue;
9506ad1575SMatthew G. Knepley for (l = 0; l < 4; ++l) {
9606ad1575SMatthew G. Knepley if ((l == i) || (l == j) || (l == k)) continue;
979566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], i, f[i]));
989566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], j, f[j]));
999566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], k, f[k]));
1009566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], l, f[l]));
1019566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 4, f, 4, f));
1029566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
10306ad1575SMatthew G. Knepley }
10406ad1575SMatthew G. Knepley }
10506ad1575SMatthew G. Knepley }
10606ad1575SMatthew G. Knepley }
10706ad1575SMatthew G. Knepley /* Check f0 and f1 in different orders */
10806ad1575SMatthew G. Knepley for (i = 0; i < 4; ++i) {
10906ad1575SMatthew G. Knepley for (j = 0; j < 4; ++j) {
11006ad1575SMatthew G. Knepley if (j == i) continue;
11106ad1575SMatthew G. Knepley for (k = 0; k < 4; ++k) {
11206ad1575SMatthew G. Knepley if ((k == i) || (k == j)) continue;
11306ad1575SMatthew G. Knepley for (l = 0; l < 4; ++l) {
11406ad1575SMatthew G. Knepley if ((l == i) || (l == j) || (l == k)) continue;
1159566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], i, f[i]));
1169566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], j, f[j]));
1179566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], k, f[k]));
1189566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], l, f[l]));
1199566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 4, f, 4, f));
1209566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
12106ad1575SMatthew G. Knepley }
12206ad1575SMatthew G. Knepley }
12306ad1575SMatthew G. Knepley }
12406ad1575SMatthew G. Knepley }
1259566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label));
1263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12706ad1575SMatthew G. Knepley }
12806ad1575SMatthew G. Knepley
TestAdd(PetscWeakForm wf)129d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestAdd(PetscWeakForm wf)
130d71ae5a4SJacob Faibussowitsch {
1312192575eSBarry Smith PetscPointFn *f[4] = {f0, f1, f2, f3}, *fp[4];
13206ad1575SMatthew G. Knepley DMLabel label;
13306ad1575SMatthew G. Knepley const PetscInt value = 3, field = 1, part = 2;
13406ad1575SMatthew G. Knepley PetscFormKey key;
13506ad1575SMatthew G. Knepley PetscInt i, j, k, l;
13606ad1575SMatthew G. Knepley
13706ad1575SMatthew G. Knepley PetscFunctionBegin;
1389566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
1399371c9d4SSatish Balay key.label = label;
1409371c9d4SSatish Balay key.value = value;
1419371c9d4SSatish Balay key.field = field;
1429371c9d4SSatish Balay key.part = part;
14306ad1575SMatthew G. Knepley /* Check f0 */
14406ad1575SMatthew G. Knepley for (i = 0; i < 4; ++i) {
14506ad1575SMatthew G. Knepley for (j = 0; j < 4; ++j) {
14606ad1575SMatthew G. Knepley if (j == i) continue;
14706ad1575SMatthew G. Knepley for (k = 0; k < 4; ++k) {
14806ad1575SMatthew G. Knepley if ((k == i) || (k == j)) continue;
14906ad1575SMatthew G. Knepley for (l = 0; l < 4; ++l) {
15006ad1575SMatthew G. Knepley if ((l == i) || (l == j) || (l == k)) continue;
1519566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[i], NULL));
1529566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[j], NULL));
1539566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[k], NULL));
1549566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[l], NULL));
1559371c9d4SSatish Balay fp[0] = f[i];
1569371c9d4SSatish Balay fp[1] = f[j];
1579371c9d4SSatish Balay fp[2] = f[k];
1589371c9d4SSatish Balay fp[3] = f[l];
1599566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 4, fp, 0, NULL));
1609566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
16106ad1575SMatthew G. Knepley }
16206ad1575SMatthew G. Knepley }
16306ad1575SMatthew G. Knepley }
16406ad1575SMatthew G. Knepley }
16506ad1575SMatthew G. Knepley /* Check f1 */
16606ad1575SMatthew G. Knepley for (i = 0; i < 4; ++i) {
16706ad1575SMatthew G. Knepley for (j = 0; j < 4; ++j) {
16806ad1575SMatthew G. Knepley if (j == i) continue;
16906ad1575SMatthew G. Knepley for (k = 0; k < 4; ++k) {
17006ad1575SMatthew G. Knepley if ((k == i) || (k == j)) continue;
17106ad1575SMatthew G. Knepley for (l = 0; l < 4; ++l) {
17206ad1575SMatthew G. Knepley if ((l == i) || (l == j) || (l == k)) continue;
1739566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[i]));
1749566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[j]));
1759566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[k]));
1769566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[l]));
1779371c9d4SSatish Balay fp[0] = f[i];
1789371c9d4SSatish Balay fp[1] = f[j];
1799371c9d4SSatish Balay fp[2] = f[k];
1809371c9d4SSatish Balay fp[3] = f[l];
1819566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 0, NULL, 4, fp));
1829566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
18306ad1575SMatthew G. Knepley }
18406ad1575SMatthew G. Knepley }
18506ad1575SMatthew G. Knepley }
18606ad1575SMatthew G. Knepley }
18706ad1575SMatthew G. Knepley /* Check f0 and f1 */
18806ad1575SMatthew G. Knepley for (i = 0; i < 4; ++i) {
18906ad1575SMatthew G. Knepley for (j = 0; j < 4; ++j) {
19006ad1575SMatthew G. Knepley if (j == i) continue;
19106ad1575SMatthew G. Knepley for (k = 0; k < 4; ++k) {
19206ad1575SMatthew G. Knepley if ((k == i) || (k == j)) continue;
19306ad1575SMatthew G. Knepley for (l = 0; l < 4; ++l) {
19406ad1575SMatthew G. Knepley if ((l == i) || (l == j) || (l == k)) continue;
1959566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[i], f[i]));
1969566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[j], f[j]));
1979566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[k], f[k]));
1989566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[l], f[l]));
1999371c9d4SSatish Balay fp[0] = f[i];
2009371c9d4SSatish Balay fp[1] = f[j];
2019371c9d4SSatish Balay fp[2] = f[k];
2029371c9d4SSatish Balay fp[3] = f[l];
2039566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 4, fp, 4, fp));
2049566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
20506ad1575SMatthew G. Knepley }
20606ad1575SMatthew G. Knepley }
20706ad1575SMatthew G. Knepley }
20806ad1575SMatthew G. Knepley }
2099566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label));
2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
21106ad1575SMatthew G. Knepley }
21206ad1575SMatthew G. Knepley
TestSetIndexAdd(PetscWeakForm wf)213d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestSetIndexAdd(PetscWeakForm wf)
214d71ae5a4SJacob Faibussowitsch {
2152192575eSBarry Smith PetscPointFn *f[4] = {f0, f1, f2, f3};
21606ad1575SMatthew G. Knepley DMLabel label;
21706ad1575SMatthew G. Knepley const PetscInt value = 3, field = 1, part = 2;
21806ad1575SMatthew G. Knepley PetscFormKey key;
21906ad1575SMatthew G. Knepley
22006ad1575SMatthew G. Knepley PetscFunctionBegin;
2219566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
2229371c9d4SSatish Balay key.label = label;
2239371c9d4SSatish Balay key.value = value;
2249371c9d4SSatish Balay key.field = field;
2259371c9d4SSatish Balay key.part = part;
22606ad1575SMatthew G. Knepley /* Check f0 */
2279566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
2289566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 1, f[1], 0, NULL));
2299566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
2309566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[3], NULL));
2319566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
2329566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
2339566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
2349566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
2359566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 2, f[2], 0, NULL));
2369566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[3], NULL));
2379566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
2389566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
2399566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
2409566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
2419566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
2429566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
2439566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
2449566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
2459566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[0], NULL));
2469566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 1, f[1], 0, NULL));
2479566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
2489566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
2499566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
2509566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
2519566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[0], NULL));
2529566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
2539566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 2, f[2], 0, NULL));
2549566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
2559566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
2569566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
25706ad1575SMatthew G. Knepley /* Check f1 */
2589566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
2599566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 1, f[1]));
2609566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
2619566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[3]));
2629566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
2639566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
2649566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
2659566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
2669566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 2, f[2]));
2679566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[3]));
2689566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
2699566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
2709566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
2719566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
2729566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
2739566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
2749566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
2759566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
2769566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[0]));
2779566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 1, f[1]));
2789566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
2799566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
2809566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
2819566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
2829566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[0]));
2839566063dSJacob Faibussowitsch PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
2849566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 2, f[2]));
2859566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
2869566063dSJacob Faibussowitsch PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
2879566063dSJacob Faibussowitsch PetscCall(PetscWeakFormClear(wf));
28806ad1575SMatthew G. Knepley
2899566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label));
2903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
29106ad1575SMatthew G. Knepley }
29206ad1575SMatthew G. Knepley
main(int argc,char ** argv)293d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
294d71ae5a4SJacob Faibussowitsch {
29506ad1575SMatthew G. Knepley PetscWeakForm wf;
29606ad1575SMatthew G. Knepley
297327415f7SBarry Smith PetscFunctionBeginUser;
2989566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2999566063dSJacob Faibussowitsch PetscCall(PetscWeakFormCreate(PETSC_COMM_SELF, &wf));
3009566063dSJacob Faibussowitsch PetscCall(TestSetIndex(wf));
3019566063dSJacob Faibussowitsch PetscCall(TestAdd(wf));
3029566063dSJacob Faibussowitsch PetscCall(TestSetIndexAdd(wf));
3039566063dSJacob Faibussowitsch PetscCall(PetscWeakFormDestroy(&wf));
3049566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
305b122ec5aSJacob Faibussowitsch return 0;
30606ad1575SMatthew G. Knepley }
30706ad1575SMatthew G. Knepley
30806ad1575SMatthew G. Knepley /*TEST
30906ad1575SMatthew G. Knepley
31006ad1575SMatthew G. Knepley test:
31106ad1575SMatthew G. Knepley suffix: 0
312*3886731fSPierre Jolivet output_file: output/empty.out
31306ad1575SMatthew G. Knepley
31406ad1575SMatthew G. Knepley TEST*/
315