xref: /petsc/include/petscdt.h (revision 8cd1e013f4f44ccfc597f967d47b0f55a2d81953)
1 /*
2   Common tools for constructing discretizations
3 */
4 #if !defined(PETSCDT_H)
5 #define PETSCDT_H
6 
7 #include <petscsys.h>
8 
9 /*S
10   PetscQuadrature - Quadrature rule for integration.
11 
12   Level: beginner
13 
14 .seealso:  PetscQuadratureCreate(), PetscQuadratureDestroy()
15 S*/
16 typedef struct _p_PetscQuadrature *PetscQuadrature;
17 
18 /*E
19   PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights
20 
21   Level: intermediate
22 
23 $  PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA - compute the nodes via linear algebra
24 $  PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method
25 
26 E*/
27 typedef enum {PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON} PetscGaussLobattoLegendreCreateType;
28 
29 PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
30 PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
31 PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*);
32 PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
33 PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*);
34 PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
35 PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]);
36 PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []);
37 PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
38 PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
39 
40 PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
41 
42 PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
43 PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
44 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*);
45 PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
46 PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
47 PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
48 
49 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
50 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
51 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
52 
53 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
54 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
55 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
56 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
57 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
58 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
59 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
60 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
61 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
62 
63 PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
64 PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
65 PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
66 PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
67 PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
68 PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
69 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
70 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
71 PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
72 
73 #if defined(PETSC_USE_64BIT_INDICES)
74 #define PETSC_FACTORIAL_MAX 20
75 #define PETSC_BINOMIAL_MAX  61
76 #else
77 #define PETSC_FACTORIAL_MAX 12
78 #define PETSC_BINOMIAL_MAX  29
79 #endif
80 
81 /*MC
82    PetscDTFactorial - Approximate n! as a real number
83 
84    Input Arguments:
85 
86 .  n - a non-negative integer
87 
88    Output Arguments;
89 
90 .  factorial - n!
91 
92    Level: beginner
93 M*/
94 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
95 {
96   PetscReal f = 1.0;
97   PetscInt  i;
98 
99   PetscFunctionBegin;
100   for (i = 1; i < n+1; ++i) f *= i;
101   *factorial = f;
102   PetscFunctionReturn(0);
103 }
104 
105 /*MC
106    PetscDTFactorialInt - Compute n! as an integer
107 
108    Input Arguments:
109 
110 .  n - a non-negative integer
111 
112    Output Arguments;
113 
114 .  factorial - n!
115 
116    Level: beginner
117 
118    Note: this is limited to n such that n! can be represented by PetscInt, which is 12 if PetscInt is a signed 32-bit integer and 20 if PetscInt is a signed 64-bit integer.
119 M*/
120 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
121 {
122   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
123 
124   PetscFunctionBeginHot;
125   if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX);
126   if (n <= 12) {
127     *factorial = facLookup[n];
128   } else {
129     PetscInt f = facLookup[12];
130     PetscInt i;
131 
132     for (i = 13; i < n+1; ++i) f *= i;
133     *factorial = f;
134   }
135   PetscFunctionReturn(0);
136 }
137 
138 /*MC
139    PetscDTBinomial - Approximate the binomial coefficient "n choose k"
140 
141    Input Arguments:
142 
143 +  n - a non-negative integer
144 -  k - an integer between 0 and n, inclusive
145 
146    Output Arguments;
147 
148 .  binomial - approximation of the binomial coefficient n choose k
149 
150    Level: beginner
151 M*/
152 PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
153 {
154   PetscFunctionBeginHot;
155   if (n < 0 || k < 0 || k > n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n\n", n, k);
156   if (n <= 3) {
157     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
158 
159     *binomial = binomLookup[n][k];
160   } else {
161     PetscReal binom = 1.;
162     PetscInt  i;
163 
164     k = PetscMin(k, n - k);
165     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
166     *binomial = binom;
167   }
168   PetscFunctionReturn(0);
169 }
170 
171 /*MC
172    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
173 
174    Input Arguments:
175 
176 +  n - a non-negative integer
177 -  k - an integer between 0 and n, inclusive
178 
179    Output Arguments;
180 
181 .  binomial - the binomial coefficient n choose k
182 
183    Note: this is limited by integers that can be represented by PetscInt
184 
185    Level: beginner
186 M*/
187 PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
188 {
189   PetscFunctionBeginHot;
190   if (n < 0 || k < 0 || k > n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n\n", n, k);
191   if (n > PETSC_BINOMIAL_MAX) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial elements %D is larger than max for PetscInt, %D\n", n, PETSC_BINOMIAL_MAX);
192   if (n <= 3) {
193     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
194 
195     *binomial = binomLookup[n][k];
196   } else {
197     PetscInt  binom = 1;
198     PetscInt  i;
199 
200     k = PetscMin(k, n - k);
201     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
202     *binomial = binom;
203   }
204   PetscFunctionReturn(0);
205 }
206 
207 /*MC
208    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
209 
210    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
211    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
212    some position j >= i.  We encode this swap as the difference (j - i).  The difference d_i at step i is less than
213    (n - i).  We encode this sequence of n-1 differences [d_0, ..., d_{n-2}] as the number
214    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
215 
216    Input Arguments:
217 
218 +  n - a non-negative integer (see note about limits below)
219 -  k - an integer in [0, n!)
220 
221    Output Arguments:
222 
223 +  perm - the permuted list of the integers [0, ..., n-1]
224 -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
225 
226    Note: this is limited to n such that n! can be represented by PetscInt, which is 12 if PetscInt is a signed 32-bit integer and 20 if PetscInt is a signed 64-bit integer.
227 
228    Level: beginner
229 M*/
230 PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
231 {
232   PetscInt  odd = 0;
233   PetscInt  i;
234   PetscInt  work[PETSC_FACTORIAL_MAX];
235   PetscInt *w;
236 
237   PetscFunctionBeginHot;
238   if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX);
239   w = &work[n - 2];
240   for (i = 2; i <= n; i++) {
241     *(w--) = k % i;
242     k /= i;
243   }
244   for (i = 0; i < n; i++) perm[i] = i;
245   for (i = 0; i < n - 1; i++) {
246     PetscInt s = work[i];
247     PetscInt swap = perm[i];
248 
249     perm[i] = perm[i + s];
250     perm[i + s] = swap;
251     odd ^= (!!s);
252   }
253   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
254   PetscFunctionReturn(0);
255 }
256 
257 /*MC
258    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts PetscDTEnumPerm.
259 
260    Input Arguments:
261 
262 +  n - a non-negative integer (see note about limits below)
263 -  perm - the permuted list of the integers [0, ..., n-1]
264 
265    Output Arguments:
266 
267 +  k - an integer in [0, n!)
268 .  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
269 
270    Note: this is limited to n such that n! can be represented by PetscInt, which is 12 if PetscInt is a signed 32-bit integer and 20 if PetscInt is a signed 64-bit integer.
271 
272    Level: beginner
273 M*/
274 PETSC_STATIC_INLINE PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
275 {
276   PetscInt  odd = 0;
277   PetscInt  i, idx;
278   PetscInt  work[PETSC_FACTORIAL_MAX];
279   PetscInt  iwork[PETSC_FACTORIAL_MAX];
280 
281   PetscFunctionBeginHot;
282   if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX);
283   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
284   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
285   for (idx = 0, i = 0; i < n - 1; i++) {
286     PetscInt j = perm[i];
287     PetscInt icur = work[i];
288     PetscInt jloc = iwork[j];
289     PetscInt diff = jloc - i;
290 
291     idx = idx * (n - i) + diff;
292     /* swap (i, jloc) */
293     work[i] = j;
294     work[jloc] = icur;
295     iwork[j] = i;
296     iwork[icur] = jloc;
297     odd ^= (!!diff);
298   }
299   *k = idx;
300   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
301   PetscFunctionReturn(0);
302 }
303 
304 /*MC
305    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
306    The encoding is in lexicographic order.
307 
308    Input Arguments:
309 
310 +  n - a non-negative integer (see note about limits below)
311 .  k - an integer in [0, n]
312 -  j - an index in [0, n choose k)
313 
314    Output Arguments:
315 
316 .  subset - the jth subset of size k of the integers [0, ..., n - 1]
317 
318    Note: this is limited by arguments such that n choose k can be represented by PetscInt
319 
320    Level: beginner
321 
322 .seealso: PetscDTSubsetIndex()
323 M*/
324 PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
325 {
326   PetscInt       Nk, i, l;
327   PetscErrorCode ierr;
328 
329   PetscFunctionBeginHot;
330   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
331   for (i = 0, l = 0; i < n && l < k; i++) {
332     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
333     PetscInt Nminusk = Nk - Nminuskminus;
334 
335     if (j < Nminuskminus) {
336       subset[l++] = i;
337       Nk = Nminuskminus;
338     } else {
339       j -= Nminuskminus;
340       Nk = Nminusk;
341     }
342   }
343   PetscFunctionReturn(0);
344 }
345 
346 /*MC
347    PetscDTSubsetIndex - Convert an ordered subset of k integers from the set [0, ..., n - 1] to its encoding as an integers in [0, n choose k) in lexicographic order.  This is the inverse of PetscDTEnumSubset.
348 
349    Input Arguments:
350 
351 +  n - a non-negative integer (see note about limits below)
352 .  k - an integer in [0, n]
353 -  subset - an ordered subset of the integers [0, ..., n - 1]
354 
355    Output Arguments:
356 
357 .  index - the rank of the subset in lexicographic order
358 
359    Note: this is limited by arguments such that n choose k can be represented by PetscInt
360 
361    Level: beginner
362 
363 .seealso: PetscDTEnumSubset()
364 M*/
365 PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
366 {
367   PetscInt       i, j = 0, l, Nk;
368   PetscErrorCode ierr;
369 
370   PetscFunctionBeginHot;
371   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
372   for (i = 0, l = 0; i < n && l < k; i++) {
373     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
374     PetscInt Nminusk = Nk - Nminuskminus;
375 
376     if (subset[l] == i) {
377       l++;
378       Nk = Nminuskminus;
379     } else {
380       j += Nminuskminus;
381       Nk = Nminusk;
382     }
383   }
384   *index = j;
385   PetscFunctionReturn(0);
386 }
387 
388 
389 /*MC
390    PetscDTEnumSubset - Split the integers [0, ..., n - 1] into two complementary ordered subsets, the first of size k and beingthe jth in lexicographic order.
391 
392    Input Arguments:
393 
394 +  n - a non-negative integer (see note about limits below)
395 .  k - an integer in [0, n]
396 -  j - an index in [0, n choose k)
397 
398    Output Arguments:
399 
400 +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
401 -  isOdd - if not NULL, return whether the permutation is even or odd.
402 
403    Note: this is limited by arguments such that n choose k can be represented by PetscInt
404 
405    Level: beginner
406 
407 .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex()
408 M*/
409 PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
410 {
411   PetscInt       i, l, m, *subcomp, Nk;
412   PetscInt       odd;
413   PetscErrorCode ierr;
414 
415   PetscFunctionBeginHot;
416   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
417   odd = 0;
418   subcomp = &perm[k];
419   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
420     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
421     PetscInt Nminusk = Nk - Nminuskminus;
422 
423     if (j < Nminuskminus) {
424       perm[l++] = i;
425       Nk = Nminuskminus;
426     } else {
427       subcomp[m++] = i;
428       j -= Nminuskminus;
429       odd ^= ((k - l) & 1);
430       Nk = Nminusk;
431     }
432   }
433   for (; i < n; i++) {
434     subcomp[m++] = i;
435   }
436   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
437   PetscFunctionReturn(0);
438 }
439 
440 #endif
441