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