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