xref: /petsc/src/dm/dt/tests/ex12.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
1 static char help[] = "Tests for PetscWeakForm.\n\n";
2 
3 #include <petscds.h>
4 
5 static void f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
6                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
7                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
8                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
9 {
10   f0[0] = 0.0;
11 }
12 
13 static void f1(PetscInt dim, PetscInt Nf, PetscInt NfAux,
14                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
15                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
16                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
17 {
18   f0[0] = 0.0;
19 }
20 
21 static void f2(PetscInt dim, PetscInt Nf, PetscInt NfAux,
22                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
23                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
24                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
25 {
26   f0[0] = 0.0;
27 }
28 
29 static void f3(PetscInt dim, PetscInt Nf, PetscInt NfAux,
30                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
31                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
32                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
33 {
34   f0[0] = 0.0;
35 }
36 
37 static PetscErrorCode CheckResidual(PetscWeakForm wf, PetscFormKey key, PetscInt in0, PetscPointFunc if0[], PetscInt in1, PetscPointFunc if1[])
38 {
39   PetscPointFunc *f0, *f1;
40   PetscInt        n0, n1, i;
41 
42   PetscFunctionBegin;
43   PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0, &n1, &f1));
44   PetscCheckFalse(n0 != in0,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Found %D f0 functions != %D functions input", n0, in0);
45   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);}
46   PetscCheckFalse(n1 != in1,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Found %D f1 functions != %D functions input", n0, in0);
47   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);}
48   PetscFunctionReturn(0);
49 }
50 
51 static PetscErrorCode TestSetIndex(PetscWeakForm wf)
52 {
53   PetscPointFunc   f[4] = {f0, f1, f2, f3};
54   DMLabel          label;
55   const PetscInt   value = 3, field = 1, part = 2;
56   PetscFormKey key;
57   PetscInt         i, j, k, l;
58 
59   PetscFunctionBegin;
60   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
61   key.label = label; key.value = value; key.field = field; key.part = part;
62   /* Check f0 */
63   for (i = 0; i < 4; ++i) {
64     for (j = 0; j < 4; ++j) {
65       if (j == i) continue;
66       for (k = 0; k < 4; ++k) {
67         if ((k == i) || (k == j)) continue;
68         for (l = 0; l < 4; ++l) {
69           if ((l == i) || (l == j) || (l == k)) continue;
70           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], 0, NULL));
71           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], 0, NULL));
72           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], 0, NULL));
73           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], 0, NULL));
74           PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
75           PetscCall(PetscWeakFormClear(wf));
76         }
77       }
78     }
79   }
80   /* Check f1 */
81   for (i = 0; i < 4; ++i) {
82     for (j = 0; j < 4; ++j) {
83       if (j == i) continue;
84       for (k = 0; k < 4; ++k) {
85         if ((k == i) || (k == j)) continue;
86         for (l = 0; l < 4; ++l) {
87           if ((l == i) || (l == j) || (l == k)) continue;
88           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, i, f[i]));
89           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, j, f[j]));
90           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, k, f[k]));
91           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, l, f[l]));
92           PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
93           PetscCall(PetscWeakFormClear(wf));
94         }
95       }
96     }
97   }
98   /* Check f0 and f1 */
99   for (i = 0; i < 4; ++i) {
100     for (j = 0; j < 4; ++j) {
101       if (j == i) continue;
102       for (k = 0; k < 4; ++k) {
103         if ((k == i) || (k == j)) continue;
104         for (l = 0; l < 4; ++l) {
105           if ((l == i) || (l == j) || (l == k)) continue;
106           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], i, f[i]));
107           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], j, f[j]));
108           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], k, f[k]));
109           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], l, f[l]));
110           PetscCall(CheckResidual(wf, key, 4, f, 4, f));
111           PetscCall(PetscWeakFormClear(wf));
112         }
113       }
114     }
115   }
116   /* Check f0 and f1 in different orders */
117   for (i = 0; i < 4; ++i) {
118     for (j = 0; j < 4; ++j) {
119       if (j == i) continue;
120       for (k = 0; k < 4; ++k) {
121         if ((k == i) || (k == j)) continue;
122         for (l = 0; l < 4; ++l) {
123           if ((l == i) || (l == j) || (l == k)) continue;
124           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, l, f[l], i, f[i]));
125           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, k, f[k], j, f[j]));
126           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, j, f[j], k, f[k]));
127           PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, i, f[i], l, f[l]));
128           PetscCall(CheckResidual(wf, key, 4, f, 4, f));
129           PetscCall(PetscWeakFormClear(wf));
130         }
131       }
132     }
133   }
134   PetscCall(DMLabelDestroy(&label));
135   PetscFunctionReturn(0);
136 }
137 
138 static PetscErrorCode TestAdd(PetscWeakForm wf)
139 {
140   PetscPointFunc   f[4] = {f0, f1, f2, f3}, fp[4];
141   DMLabel          label;
142   const PetscInt   value = 3, field = 1, part = 2;
143   PetscFormKey key;
144   PetscInt         i, j, k, l;
145 
146   PetscFunctionBegin;
147   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
148   key.label = label; key.value = value; key.field = field; key.part = part;
149   /* Check f0 */
150   for (i = 0; i < 4; ++i) {
151     for (j = 0; j < 4; ++j) {
152       if (j == i) continue;
153       for (k = 0; k < 4; ++k) {
154         if ((k == i) || (k == j)) continue;
155         for (l = 0; l < 4; ++l) {
156           if ((l == i) || (l == j) || (l == k)) continue;
157           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[i], NULL));
158           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[j], NULL));
159           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[k], NULL));
160           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[l], NULL));
161           fp[0] = f[i]; fp[1] = f[j]; fp[2] = f[k]; fp[3] = f[l];
162           PetscCall(CheckResidual(wf, key, 4, fp, 0, NULL));
163           PetscCall(PetscWeakFormClear(wf));
164         }
165       }
166     }
167   }
168   /* Check f1 */
169   for (i = 0; i < 4; ++i) {
170     for (j = 0; j < 4; ++j) {
171       if (j == i) continue;
172       for (k = 0; k < 4; ++k) {
173         if ((k == i) || (k == j)) continue;
174         for (l = 0; l < 4; ++l) {
175           if ((l == i) || (l == j) || (l == k)) continue;
176           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[i]));
177           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[j]));
178           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[k]));
179           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[l]));
180           fp[0] = f[i]; fp[1] = f[j]; fp[2] = f[k]; fp[3] = f[l];
181           PetscCall(CheckResidual(wf, key, 0, NULL, 4, fp));
182           PetscCall(PetscWeakFormClear(wf));
183         }
184       }
185     }
186   }
187   /* Check f0 and f1 */
188   for (i = 0; i < 4; ++i) {
189     for (j = 0; j < 4; ++j) {
190       if (j == i) continue;
191       for (k = 0; k < 4; ++k) {
192         if ((k == i) || (k == j)) continue;
193         for (l = 0; l < 4; ++l) {
194           if ((l == i) || (l == j) || (l == k)) continue;
195           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[i], f[i]));
196           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[j], f[j]));
197           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[k], f[k]));
198           PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[l], f[l]));
199           fp[0] = f[i]; fp[1] = f[j]; fp[2] = f[k]; fp[3] = f[l];
200           PetscCall(CheckResidual(wf, key, 4, fp, 4, fp));
201           PetscCall(PetscWeakFormClear(wf));
202         }
203       }
204     }
205   }
206   PetscCall(DMLabelDestroy(&label));
207   PetscFunctionReturn(0);
208 }
209 
210 static PetscErrorCode TestSetIndexAdd(PetscWeakForm wf)
211 {
212   PetscPointFunc   f[4] = {f0, f1, f2, f3};
213   DMLabel          label;
214   const PetscInt   value = 3, field = 1, part = 2;
215   PetscFormKey key;
216 
217   PetscFunctionBegin;
218   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Test", &label));
219   key.label = label; key.value = value; key.field = field; key.part = part;
220   /* Check f0 */
221   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
222   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 1, f[1], 0, NULL));
223   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
224   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[3], NULL));
225   PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
226   PetscCall(PetscWeakFormClear(wf));
227   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
228   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
229   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 2, f[2], 0, NULL));
230   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[3], NULL));
231   PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
232   PetscCall(PetscWeakFormClear(wf));
233   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, f[0], 0, NULL));
234   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
235   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
236   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
237   PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
238   PetscCall(PetscWeakFormClear(wf));
239   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[0], NULL));
240   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 1, f[1], 0, NULL));
241   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[2], NULL));
242   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
243   PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
244   PetscCall(PetscWeakFormClear(wf));
245   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[0], NULL));
246   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, f[1], NULL));
247   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 2, f[2], 0, NULL));
248   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 3, f[3], 0, NULL));
249   PetscCall(CheckResidual(wf, key, 4, f, 0, NULL));
250   PetscCall(PetscWeakFormClear(wf));
251   /* Check f1 */
252   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
253   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 1, f[1]));
254   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
255   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[3]));
256   PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
257   PetscCall(PetscWeakFormClear(wf));
258   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
259   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
260   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 2, f[2]));
261   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[3]));
262   PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
263   PetscCall(PetscWeakFormClear(wf));
264   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 0, f[0]));
265   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
266   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
267   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
268   PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
269   PetscCall(PetscWeakFormClear(wf));
270   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[0]));
271   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 1, f[1]));
272   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[2]));
273   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
274   PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
275   PetscCall(PetscWeakFormClear(wf));
276   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[0]));
277   PetscCall(PetscWeakFormAddResidual(wf, key.label, key.value, key.field, key.part, NULL, f[1]));
278   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 2, f[2]));
279   PetscCall(PetscWeakFormSetIndexResidual(wf, key.label, key.value, key.field, key.part, 0, NULL, 3, f[3]));
280   PetscCall(CheckResidual(wf, key, 0, NULL, 4, f));
281   PetscCall(PetscWeakFormClear(wf));
282 
283   PetscCall(DMLabelDestroy(&label));
284   PetscFunctionReturn(0);
285 }
286 
287 int main(int argc,char **argv)
288 {
289   PetscWeakForm  wf;
290 
291   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
292   PetscCall(PetscWeakFormCreate(PETSC_COMM_SELF, &wf));
293   PetscCall(TestSetIndex(wf));
294   PetscCall(TestAdd(wf));
295   PetscCall(TestSetIndexAdd(wf));
296   PetscCall(PetscWeakFormDestroy(&wf));
297   PetscCall(PetscFinalize());
298   return 0;
299 }
300 
301 /*TEST
302 
303   test:
304     suffix: 0
305 
306 TEST*/
307