xref: /petsc/src/dm/dt/tests/ex12.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
106ad1575SMatthew G. Knepley static char help[] = "Tests for PetscWeakForm.\n\n";
206ad1575SMatthew G. Knepley 
306ad1575SMatthew G. Knepley #include <petscds.h>
406ad1575SMatthew G. Knepley 
506ad1575SMatthew G. Knepley static void f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
606ad1575SMatthew G. Knepley                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
706ad1575SMatthew G. Knepley                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
806ad1575SMatthew G. Knepley                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
906ad1575SMatthew G. Knepley {
1006ad1575SMatthew G. Knepley   f0[0] = 0.0;
1106ad1575SMatthew G. Knepley }
1206ad1575SMatthew G. Knepley 
1306ad1575SMatthew G. Knepley static void f1(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1406ad1575SMatthew G. Knepley                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1506ad1575SMatthew G. Knepley                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1606ad1575SMatthew G. Knepley                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1706ad1575SMatthew G. Knepley {
1806ad1575SMatthew G. Knepley   f0[0] = 0.0;
1906ad1575SMatthew G. Knepley }
2006ad1575SMatthew G. Knepley 
2106ad1575SMatthew G. Knepley static void f2(PetscInt dim, PetscInt Nf, PetscInt NfAux,
2206ad1575SMatthew G. Knepley                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
2306ad1575SMatthew G. Knepley                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
2406ad1575SMatthew G. Knepley                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
2506ad1575SMatthew G. Knepley {
2606ad1575SMatthew G. Knepley   f0[0] = 0.0;
2706ad1575SMatthew G. Knepley }
2806ad1575SMatthew G. Knepley 
2906ad1575SMatthew G. Knepley static void f3(PetscInt dim, PetscInt Nf, PetscInt NfAux,
3006ad1575SMatthew G. Knepley                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
3106ad1575SMatthew G. Knepley                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
3206ad1575SMatthew G. Knepley                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
3306ad1575SMatthew G. Knepley {
3406ad1575SMatthew G. Knepley   f0[0] = 0.0;
3506ad1575SMatthew G. Knepley }
3606ad1575SMatthew G. Knepley 
3706ad1575SMatthew G. Knepley static PetscErrorCode CheckResidual(PetscWeakForm wf, PetscFormKey key, PetscInt in0, PetscPointFunc if0[], PetscInt in1, PetscPointFunc if1[])
3806ad1575SMatthew G. Knepley {
3906ad1575SMatthew G. Knepley   PetscPointFunc *f0, *f1;
4006ad1575SMatthew G. Knepley   PetscInt        n0, n1, i;
4106ad1575SMatthew G. Knepley 
4206ad1575SMatthew G. Knepley   PetscFunctionBegin;
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0, &n1, &f1));
442c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n0 != in0,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Found %D f0 functions != %D functions input", n0, in0);
452c71b3e2SJacob Faibussowitsch   for (i = 0; i < n0; ++i) {PetscCheckFalse(f0[i] != if0[i],PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "f0[%D] != input f0[%D]", i, i);}
462c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n1 != in1,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Found %D f1 functions != %D functions input", n0, in0);
472c71b3e2SJacob Faibussowitsch   for (i = 0; i < n1; ++i) {PetscCheckFalse(f1[i] != if1[i],PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "f1[%D] != input f1[%D]", i, i);}
4806ad1575SMatthew G. Knepley   PetscFunctionReturn(0);
4906ad1575SMatthew G. Knepley }
5006ad1575SMatthew G. Knepley 
5106ad1575SMatthew G. Knepley static PetscErrorCode TestSetIndex(PetscWeakForm wf)
5206ad1575SMatthew G. Knepley {
5306ad1575SMatthew G. Knepley   PetscPointFunc   f[4] = {f0, f1, f2, f3};
5406ad1575SMatthew G. Knepley   DMLabel          label;
5506ad1575SMatthew G. Knepley   const PetscInt   value = 3, field = 1, part = 2;
5606ad1575SMatthew G. Knepley   PetscFormKey key;
5706ad1575SMatthew G. Knepley   PetscInt         i, j, k, l;
5806ad1575SMatthew G. Knepley 
5906ad1575SMatthew G. Knepley   PetscFunctionBegin;
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
6106ad1575SMatthew G. Knepley   key.label = label; key.value = value; key.field = field; key.part = part;
6206ad1575SMatthew G. Knepley   /* Check f0 */
6306ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
6406ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
6506ad1575SMatthew G. Knepley       if (j == i) continue;
6606ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
6706ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
6806ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
6906ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
70*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], 0, NULL));
71*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], 0, NULL));
72*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], 0, NULL));
73*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], 0, NULL));
74*5f80ce2aSJacob Faibussowitsch           CHKERRQ(CheckResidual(wf, key, 4, f, 0, NULL));
75*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormClear(wf));
7606ad1575SMatthew G. Knepley         }
7706ad1575SMatthew G. Knepley       }
7806ad1575SMatthew G. Knepley     }
7906ad1575SMatthew G. Knepley   }
8006ad1575SMatthew G. Knepley   /* Check f1 */
8106ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
8206ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
8306ad1575SMatthew G. Knepley       if (j == i) continue;
8406ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
8506ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
8606ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
8706ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
88*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, i, f[i]));
89*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, j, f[j]));
90*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, k, f[k]));
91*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, l, f[l]));
92*5f80ce2aSJacob Faibussowitsch           CHKERRQ(CheckResidual(wf, key, 0, NULL, 4, f));
93*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormClear(wf));
9406ad1575SMatthew G. Knepley         }
9506ad1575SMatthew G. Knepley       }
9606ad1575SMatthew G. Knepley     }
9706ad1575SMatthew G. Knepley   }
9806ad1575SMatthew G. Knepley   /* Check f0 and f1 */
9906ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
10006ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
10106ad1575SMatthew G. Knepley       if (j == i) continue;
10206ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
10306ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
10406ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
10506ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
106*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], i, f[i]));
107*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], j, f[j]));
108*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], k, f[k]));
109*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], l, f[l]));
110*5f80ce2aSJacob Faibussowitsch           CHKERRQ(CheckResidual(wf, key, 4, f, 4, f));
111*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormClear(wf));
11206ad1575SMatthew G. Knepley         }
11306ad1575SMatthew G. Knepley       }
11406ad1575SMatthew G. Knepley     }
11506ad1575SMatthew G. Knepley   }
11606ad1575SMatthew G. Knepley   /* Check f0 and f1 in different orders */
11706ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
11806ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
11906ad1575SMatthew G. Knepley       if (j == i) continue;
12006ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
12106ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
12206ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
12306ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
124*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], i, f[i]));
125*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], j, f[j]));
126*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], k, f[k]));
127*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], l, f[l]));
128*5f80ce2aSJacob Faibussowitsch           CHKERRQ(CheckResidual(wf, key, 4, f, 4, f));
129*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormClear(wf));
13006ad1575SMatthew G. Knepley         }
13106ad1575SMatthew G. Knepley       }
13206ad1575SMatthew G. Knepley     }
13306ad1575SMatthew G. Knepley   }
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelDestroy(&label));
13506ad1575SMatthew G. Knepley   PetscFunctionReturn(0);
13606ad1575SMatthew G. Knepley }
13706ad1575SMatthew G. Knepley 
13806ad1575SMatthew G. Knepley static PetscErrorCode TestAdd(PetscWeakForm wf)
13906ad1575SMatthew G. Knepley {
14006ad1575SMatthew G. Knepley   PetscPointFunc   f[4] = {f0, f1, f2, f3}, fp[4];
14106ad1575SMatthew G. Knepley   DMLabel          label;
14206ad1575SMatthew G. Knepley   const PetscInt   value = 3, field = 1, part = 2;
14306ad1575SMatthew G. Knepley   PetscFormKey key;
14406ad1575SMatthew G. Knepley   PetscInt         i, j, k, l;
14506ad1575SMatthew G. Knepley 
14606ad1575SMatthew G. Knepley   PetscFunctionBegin;
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
14806ad1575SMatthew G. Knepley   key.label = label; key.value = value; key.field = field; key.part = part;
14906ad1575SMatthew G. Knepley   /* Check f0 */
15006ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
15106ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
15206ad1575SMatthew G. Knepley       if (j == i) continue;
15306ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
15406ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
15506ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
15606ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
157*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[i], NULL));
158*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[j], NULL));
159*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[k], NULL));
160*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[l], NULL));
16106ad1575SMatthew G. Knepley           fp[0] = f[i]; fp[1] = f[j]; fp[2] = f[k]; fp[3] = f[l];
162*5f80ce2aSJacob Faibussowitsch           CHKERRQ(CheckResidual(wf, key, 4, fp, 0, NULL));
163*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormClear(wf));
16406ad1575SMatthew G. Knepley         }
16506ad1575SMatthew G. Knepley       }
16606ad1575SMatthew G. Knepley     }
16706ad1575SMatthew G. Knepley   }
16806ad1575SMatthew G. Knepley   /* Check f1 */
16906ad1575SMatthew G. Knepley   for (i = 0; i < 4; ++i) {
17006ad1575SMatthew G. Knepley     for (j = 0; j < 4; ++j) {
17106ad1575SMatthew G. Knepley       if (j == i) continue;
17206ad1575SMatthew G. Knepley       for (k = 0; k < 4; ++k) {
17306ad1575SMatthew G. Knepley         if ((k == i) || (k == j)) continue;
17406ad1575SMatthew G. Knepley         for (l = 0; l < 4; ++l) {
17506ad1575SMatthew G. Knepley           if ((l == i) || (l == j) || (l == k)) continue;
176*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[i]));
177*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[j]));
178*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[k]));
179*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[l]));
18006ad1575SMatthew G. Knepley           fp[0] = f[i]; fp[1] = f[j]; fp[2] = f[k]; fp[3] = f[l];
181*5f80ce2aSJacob Faibussowitsch           CHKERRQ(CheckResidual(wf, key, 0, NULL, 4, fp));
182*5f80ce2aSJacob Faibussowitsch           CHKERRQ(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;
195*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[i], f[i]));
196*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[j], f[j]));
197*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[k], f[k]));
198*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[l], f[l]));
19906ad1575SMatthew G. Knepley           fp[0] = f[i]; fp[1] = f[j]; fp[2] = f[k]; fp[3] = f[l];
200*5f80ce2aSJacob Faibussowitsch           CHKERRQ(CheckResidual(wf, key, 4, fp, 4, fp));
201*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscWeakFormClear(wf));
20206ad1575SMatthew G. Knepley         }
20306ad1575SMatthew G. Knepley       }
20406ad1575SMatthew G. Knepley     }
20506ad1575SMatthew G. Knepley   }
206*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelDestroy(&label));
20706ad1575SMatthew G. Knepley   PetscFunctionReturn(0);
20806ad1575SMatthew G. Knepley }
20906ad1575SMatthew G. Knepley 
21006ad1575SMatthew G. Knepley static PetscErrorCode TestSetIndexAdd(PetscWeakForm wf)
21106ad1575SMatthew G. Knepley {
21206ad1575SMatthew G. Knepley   PetscPointFunc   f[4] = {f0, f1, f2, f3};
21306ad1575SMatthew G. Knepley   DMLabel          label;
21406ad1575SMatthew G. Knepley   const PetscInt   value = 3, field = 1, part = 2;
21506ad1575SMatthew G. Knepley   PetscFormKey key;
21606ad1575SMatthew G. Knepley 
21706ad1575SMatthew G. Knepley   PetscFunctionBegin;
218*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
21906ad1575SMatthew G. Knepley   key.label = label; key.value = value; key.field = field; key.part = part;
22006ad1575SMatthew G. Knepley   /* Check f0 */
221*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 1, f[1], 0, NULL));
223*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
224*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[3], NULL));
225*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckResidual(wf, key, 4, f, 0, NULL));
226*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormClear(wf));
227*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
228*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
229*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 2, f[2], 0, NULL));
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[3], NULL));
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckResidual(wf, key, 4, f, 0, NULL));
232*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormClear(wf));
233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
234*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
235*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
236*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
237*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckResidual(wf, key, 4, f, 0, NULL));
238*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormClear(wf));
239*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[0], NULL));
240*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 1, f[1], 0, NULL));
241*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
242*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
243*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckResidual(wf, key, 4, f, 0, NULL));
244*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormClear(wf));
245*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[0], NULL));
246*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
247*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 2, f[2], 0, NULL));
248*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
249*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckResidual(wf, key, 4, f, 0, NULL));
250*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormClear(wf));
25106ad1575SMatthew G. Knepley   /* Check f1 */
252*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
253*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 1, f[1]));
254*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
255*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[3]));
256*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckResidual(wf, key, 0, NULL, 4, f));
257*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormClear(wf));
258*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
259*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
260*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 2, f[2]));
261*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[3]));
262*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckResidual(wf, key, 0, NULL, 4, f));
263*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormClear(wf));
264*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
265*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
266*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
267*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
268*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckResidual(wf, key, 0, NULL, 4, f));
269*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormClear(wf));
270*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[0]));
271*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 1, f[1]));
272*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
273*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
274*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckResidual(wf, key, 0, NULL, 4, f));
275*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormClear(wf));
276*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[0]));
277*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
278*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 2, f[2]));
279*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
280*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckResidual(wf, key, 0, NULL, 4, f));
281*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormClear(wf));
28206ad1575SMatthew G. Knepley 
283*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelDestroy(&label));
28406ad1575SMatthew G. Knepley   PetscFunctionReturn(0);
28506ad1575SMatthew G. Knepley }
28606ad1575SMatthew G. Knepley 
28706ad1575SMatthew G. Knepley int main(int argc,char **argv)
28806ad1575SMatthew G. Knepley {
28906ad1575SMatthew G. Knepley   PetscWeakForm  wf;
29006ad1575SMatthew G. Knepley   PetscErrorCode ierr;
29106ad1575SMatthew G. Knepley 
29206ad1575SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
293*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormCreate(PETSC_COMM_SELF, &wf));
294*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TestSetIndex(wf));
295*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TestAdd(wf));
296*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TestSetIndexAdd(wf));
297*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormDestroy(&wf));
29806ad1575SMatthew G. Knepley   ierr = PetscFinalize();
29906ad1575SMatthew G. Knepley   return ierr;
30006ad1575SMatthew G. Knepley }
30106ad1575SMatthew G. Knepley 
30206ad1575SMatthew G. Knepley /*TEST
30306ad1575SMatthew G. Knepley 
30406ad1575SMatthew G. Knepley   test:
30506ad1575SMatthew G. Knepley     suffix: 0
30606ad1575SMatthew G. Knepley 
30706ad1575SMatthew G. Knepley TEST*/
308