xref: /petsc/src/dm/dt/tests/ex12.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
106ad1575SMatthew G. Knepley static char help[] = "Tests for PetscWeakForm.\n\n";
206ad1575SMatthew G. Knepley 
306ad1575SMatthew G. Knepley #include <petscds.h>
406ad1575SMatthew G. Knepley 
5*9371c9d4SSatish Balay 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[]) {
606ad1575SMatthew G. Knepley   f0[0] = 0.0;
706ad1575SMatthew G. Knepley }
806ad1575SMatthew G. Knepley 
9*9371c9d4SSatish Balay 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[]) {
1006ad1575SMatthew G. Knepley   f0[0] = 0.0;
1106ad1575SMatthew G. Knepley }
1206ad1575SMatthew G. Knepley 
13*9371c9d4SSatish Balay 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[]) {
1406ad1575SMatthew G. Knepley   f0[0] = 0.0;
1506ad1575SMatthew G. Knepley }
1606ad1575SMatthew G. Knepley 
17*9371c9d4SSatish Balay 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[]) {
1806ad1575SMatthew G. Knepley   f0[0] = 0.0;
1906ad1575SMatthew G. Knepley }
2006ad1575SMatthew G. Knepley 
21*9371c9d4SSatish Balay static PetscErrorCode CheckResidual(PetscWeakForm wf, PetscFormKey key, PetscInt in0, PetscPointFunc if0[], PetscInt in1, PetscPointFunc if1[]) {
2206ad1575SMatthew G. Knepley   PetscPointFunc *f0, *f1;
2306ad1575SMatthew G. Knepley   PetscInt        n0, n1, i;
2406ad1575SMatthew G. Knepley 
2506ad1575SMatthew G. Knepley   PetscFunctionBegin;
269566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0, &n1, &f1));
2763a3b9bcSJacob Faibussowitsch   PetscCheck(n0 == in0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Found %" PetscInt_FMT " f0 functions != %" PetscInt_FMT " functions input", n0, in0);
2863a3b9bcSJacob Faibussowitsch   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); }
2963a3b9bcSJacob Faibussowitsch   PetscCheck(n1 == in1, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Found %" PetscInt_FMT " f1 functions != %" PetscInt_FMT " functions input", n0, in0);
3063a3b9bcSJacob Faibussowitsch   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); }
3106ad1575SMatthew G. Knepley   PetscFunctionReturn(0);
3206ad1575SMatthew G. Knepley }
3306ad1575SMatthew G. Knepley 
34*9371c9d4SSatish Balay static PetscErrorCode TestSetIndex(PetscWeakForm wf) {
3506ad1575SMatthew G. Knepley   PetscPointFunc f[4] = {f0, f1, f2, f3};
3606ad1575SMatthew G. Knepley   DMLabel        label;
3706ad1575SMatthew G. Knepley   const PetscInt value = 3, field = 1, part = 2;
3806ad1575SMatthew G. Knepley   PetscFormKey   key;
3906ad1575SMatthew G. Knepley   PetscInt       i, j, k, l;
4006ad1575SMatthew G. Knepley 
4106ad1575SMatthew G. Knepley   PetscFunctionBegin;
429566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
43*9371c9d4SSatish Balay   key.label = label;
44*9371c9d4SSatish Balay   key.value = value;
45*9371c9d4SSatish Balay   key.field = field;
46*9371c9d4SSatish Balay   key.part  = part;
4706ad1575SMatthew G. Knepley   /* Check f0 */
4806ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
4906ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
5006ad1575SMatthew G. Knepley       if (j == i) continue;
5106ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
5206ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
5306ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
5406ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
559566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], 0, NULL));
569566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], 0, NULL));
579566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], 0, NULL));
589566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], 0, NULL));
599566063dSJacob Faibussowitsch           PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
609566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormClear(wf));
6106ad1575SMatthew G. Knepley         }
6206ad1575SMatthew G. Knepley       }
6306ad1575SMatthew G. Knepley     }
6406ad1575SMatthew G. Knepley   }
6506ad1575SMatthew G. Knepley   /* Check f1 */
6606ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
6706ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
6806ad1575SMatthew G. Knepley       if (j == i) continue;
6906ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
7006ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
7106ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
7206ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
739566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, i, f[i]));
749566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, j, f[j]));
759566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, k, f[k]));
769566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, l, f[l]));
779566063dSJacob Faibussowitsch           PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
789566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormClear(wf));
7906ad1575SMatthew G. Knepley         }
8006ad1575SMatthew G. Knepley       }
8106ad1575SMatthew G. Knepley     }
8206ad1575SMatthew G. Knepley   }
8306ad1575SMatthew G. Knepley   /* Check f0 and f1 */
8406ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
8506ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
8606ad1575SMatthew G. Knepley       if (j == i) continue;
8706ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
8806ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
8906ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
9006ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
919566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], i, f[i]));
929566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], j, f[j]));
939566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], k, f[k]));
949566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], l, f[l]));
959566063dSJacob Faibussowitsch           PetscCall(CheckResidual(wf, key, 4, f, 4, f));
969566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormClear(wf));
9706ad1575SMatthew G. Knepley         }
9806ad1575SMatthew G. Knepley       }
9906ad1575SMatthew G. Knepley     }
10006ad1575SMatthew G. Knepley   }
10106ad1575SMatthew G. Knepley   /* Check f0 and f1 in different orders */
10206ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
10306ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
10406ad1575SMatthew G. Knepley       if (j == i) continue;
10506ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
10606ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
10706ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
10806ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
1099566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], i, f[i]));
1109566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], j, f[j]));
1119566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], k, f[k]));
1129566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], l, f[l]));
1139566063dSJacob Faibussowitsch           PetscCall(CheckResidual(wf, key, 4, f, 4, f));
1149566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormClear(wf));
11506ad1575SMatthew G. Knepley         }
11606ad1575SMatthew G. Knepley       }
11706ad1575SMatthew G. Knepley     }
11806ad1575SMatthew G. Knepley   }
1199566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
12006ad1575SMatthew G. Knepley   PetscFunctionReturn(0);
12106ad1575SMatthew G. Knepley }
12206ad1575SMatthew G. Knepley 
123*9371c9d4SSatish Balay static PetscErrorCode TestAdd(PetscWeakForm wf) {
12406ad1575SMatthew G. Knepley   PetscPointFunc f[4] = {f0, f1, f2, f3}, fp[4];
12506ad1575SMatthew G. Knepley   DMLabel        label;
12606ad1575SMatthew G. Knepley   const PetscInt value = 3, field = 1, part = 2;
12706ad1575SMatthew G. Knepley   PetscFormKey   key;
12806ad1575SMatthew G. Knepley   PetscInt       i, j, k, l;
12906ad1575SMatthew G. Knepley 
13006ad1575SMatthew G. Knepley   PetscFunctionBegin;
1319566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
132*9371c9d4SSatish Balay   key.label = label;
133*9371c9d4SSatish Balay   key.value = value;
134*9371c9d4SSatish Balay   key.field = field;
135*9371c9d4SSatish Balay   key.part  = part;
13606ad1575SMatthew G. Knepley   /* Check f0 */
13706ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
13806ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
13906ad1575SMatthew G. Knepley       if (j == i) continue;
14006ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
14106ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
14206ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
14306ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
1449566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[i], NULL));
1459566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[j], NULL));
1469566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[k], NULL));
1479566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[l], NULL));
148*9371c9d4SSatish Balay           fp[0] = f[i];
149*9371c9d4SSatish Balay           fp[1] = f[j];
150*9371c9d4SSatish Balay           fp[2] = f[k];
151*9371c9d4SSatish Balay           fp[3] = f[l];
1529566063dSJacob Faibussowitsch           PetscCall(CheckResidual(wf, key, 4, fp, 0, NULL));
1539566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormClear(wf));
15406ad1575SMatthew G. Knepley         }
15506ad1575SMatthew G. Knepley       }
15606ad1575SMatthew G. Knepley     }
15706ad1575SMatthew G. Knepley   }
15806ad1575SMatthew G. Knepley   /* Check f1 */
15906ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
16006ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
16106ad1575SMatthew G. Knepley       if (j == i) continue;
16206ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
16306ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
16406ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
16506ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
1669566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[i]));
1679566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[j]));
1689566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[k]));
1699566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[l]));
170*9371c9d4SSatish Balay           fp[0] = f[i];
171*9371c9d4SSatish Balay           fp[1] = f[j];
172*9371c9d4SSatish Balay           fp[2] = f[k];
173*9371c9d4SSatish Balay           fp[3] = f[l];
1749566063dSJacob Faibussowitsch           PetscCall(CheckResidual(wf, key, 0, NULL, 4, fp));
1759566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormClear(wf));
17606ad1575SMatthew G. Knepley         }
17706ad1575SMatthew G. Knepley       }
17806ad1575SMatthew G. Knepley     }
17906ad1575SMatthew G. Knepley   }
18006ad1575SMatthew G. Knepley   /* Check f0 and f1 */
18106ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
18206ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
18306ad1575SMatthew G. Knepley       if (j == i) continue;
18406ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
18506ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
18606ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
18706ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
1889566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[i], f[i]));
1899566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[j], f[j]));
1909566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[k], f[k]));
1919566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[l], f[l]));
192*9371c9d4SSatish Balay           fp[0] = f[i];
193*9371c9d4SSatish Balay           fp[1] = f[j];
194*9371c9d4SSatish Balay           fp[2] = f[k];
195*9371c9d4SSatish Balay           fp[3] = f[l];
1969566063dSJacob Faibussowitsch           PetscCall(CheckResidual(wf, key, 4, fp, 4, fp));
1979566063dSJacob Faibussowitsch           PetscCall(PetscWeakFormClear(wf));
19806ad1575SMatthew G. Knepley         }
19906ad1575SMatthew G. Knepley       }
20006ad1575SMatthew G. Knepley     }
20106ad1575SMatthew G. Knepley   }
2029566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
20306ad1575SMatthew G. Knepley   PetscFunctionReturn(0);
20406ad1575SMatthew G. Knepley }
20506ad1575SMatthew G. Knepley 
206*9371c9d4SSatish Balay static PetscErrorCode TestSetIndexAdd(PetscWeakForm wf) {
20706ad1575SMatthew G. Knepley   PetscPointFunc f[4] = {f0, f1, f2, f3};
20806ad1575SMatthew G. Knepley   DMLabel        label;
20906ad1575SMatthew G. Knepley   const PetscInt value = 3, field = 1, part = 2;
21006ad1575SMatthew G. Knepley   PetscFormKey   key;
21106ad1575SMatthew G. Knepley 
21206ad1575SMatthew G. Knepley   PetscFunctionBegin;
2139566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
214*9371c9d4SSatish Balay   key.label = label;
215*9371c9d4SSatish Balay   key.value = value;
216*9371c9d4SSatish Balay   key.field = field;
217*9371c9d4SSatish Balay   key.part  = part;
21806ad1575SMatthew G. Knepley   /* Check f0 */
2199566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
2209566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 1, f[1], 0, NULL));
2219566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
2229566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[3], NULL));
2239566063dSJacob Faibussowitsch   PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
2249566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormClear(wf));
2259566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
2269566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
2279566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 2, f[2], 0, NULL));
2289566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[3], NULL));
2299566063dSJacob Faibussowitsch   PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
2309566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormClear(wf));
2319566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
2329566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
2339566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
2349566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
2359566063dSJacob Faibussowitsch   PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
2369566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormClear(wf));
2379566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[0], NULL));
2389566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 1, f[1], 0, NULL));
2399566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
2409566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
2419566063dSJacob Faibussowitsch   PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
2429566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormClear(wf));
2439566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[0], NULL));
2449566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
2459566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 2, f[2], 0, NULL));
2469566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
2479566063dSJacob Faibussowitsch   PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
2489566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormClear(wf));
24906ad1575SMatthew G. Knepley   /* Check f1 */
2509566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
2519566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 1, f[1]));
2529566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
2539566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[3]));
2549566063dSJacob Faibussowitsch   PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
2559566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormClear(wf));
2569566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
2579566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
2589566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 2, f[2]));
2599566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[3]));
2609566063dSJacob Faibussowitsch   PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
2619566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormClear(wf));
2629566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
2639566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
2649566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
2659566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
2669566063dSJacob Faibussowitsch   PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
2679566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormClear(wf));
2689566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[0]));
2699566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 1, f[1]));
2709566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
2719566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
2729566063dSJacob Faibussowitsch   PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
2739566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormClear(wf));
2749566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[0]));
2759566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
2769566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 2, f[2]));
2779566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
2789566063dSJacob Faibussowitsch   PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
2799566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormClear(wf));
28006ad1575SMatthew G. Knepley 
2819566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
28206ad1575SMatthew G. Knepley   PetscFunctionReturn(0);
28306ad1575SMatthew G. Knepley }
28406ad1575SMatthew G. Knepley 
285*9371c9d4SSatish Balay int main(int argc, char **argv) {
28606ad1575SMatthew G. Knepley   PetscWeakForm wf;
28706ad1575SMatthew G. Knepley 
288327415f7SBarry Smith   PetscFunctionBeginUser;
2899566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2909566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormCreate(PETSC_COMM_SELF, &wf));
2919566063dSJacob Faibussowitsch   PetscCall(TestSetIndex(wf));
2929566063dSJacob Faibussowitsch   PetscCall(TestAdd(wf));
2939566063dSJacob Faibussowitsch   PetscCall(TestSetIndexAdd(wf));
2949566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormDestroy(&wf));
2959566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
296b122ec5aSJacob Faibussowitsch   return 0;
29706ad1575SMatthew G. Knepley }
29806ad1575SMatthew G. Knepley 
29906ad1575SMatthew G. Knepley /*TEST
30006ad1575SMatthew G. Knepley 
30106ad1575SMatthew G. Knepley   test:
30206ad1575SMatthew G. Knepley     suffix: 0
30306ad1575SMatthew G. Knepley 
30406ad1575SMatthew G. Knepley TEST*/
305