Skip to content

Commit 83bf834

Browse files
committed
Hypergeometric2F1
1 parent f65a287 commit 83bf834

4 files changed

Lines changed: 50 additions & 0 deletions

File tree

ufl/algorithms/apply_derivatives.py

Lines changed: 14 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -90,6 +90,7 @@
9090
Cosh,
9191
Erf,
9292
Exp,
93+
Hypergeometric2F1,
9394
Ln,
9495
MathFunction,
9596
Sin,
@@ -583,6 +584,19 @@ def _(self, o: Expr, fp: Expr) -> Expr:
583584
(f,) = o.ufl_operands
584585
return fp * (2.0 / sqrt(pi) * exp(-(f**2)))
585586

587+
@process.register(Hypergeometric2F1)
588+
@DAGTraverser.postorder
589+
def _(self, o: Expr, ap: Expr, bp: Expr, cp: Expr, fp: Expr) -> Expr:
590+
"""Differentiate a bessel_j."""
591+
a, b, c, f = o.ufl_operands
592+
if not all(nup is None or isinstance(nup, Zero) for nup in (ap, bp, cp)):
593+
raise NotImplementedError(
594+
"Differentiation of bessel function w.r.t. (a, b, c) is not supported."
595+
)
596+
597+
op = (a*b/c) * Hypergeometric2F1(a+1, b+1, c+1, f)
598+
return op * fp
599+
586600
# --- Bessel functions
587601

588602
@process.register(BesselJ)

ufl/classes.py

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -166,6 +166,7 @@
166166
Cosh,
167167
Erf,
168168
Exp,
169+
Hypergeometric2F1,
169170
Ln,
170171
MathFunction,
171172
Sin,

ufl/mathfunctions.py

Lines changed: 34 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -388,6 +388,40 @@ def evaluate(self, x, mapping, component, index_values):
388388
return math.erf(a)
389389

390390

391+
@ufl_type(is_scalar=True, num_ops=4)
392+
class Hypergeometric2F1(Operator):
393+
"""Class for hypergeometric 2F1 functions."""
394+
395+
__slots__ = "_name"
396+
397+
def __init__(self, a, b, c, argument):
398+
"""Initialise."""
399+
if not all(is_true_ufl_scalar(f) for f in (a, b, c)):
400+
raise ValueError("Expecting scalar (a, b, c).")
401+
if not is_true_ufl_scalar(argument):
402+
raise ValueError("Expecting scalar argument.")
403+
404+
Operator.__init__(self, (a, b, c, argument))
405+
406+
self._name = "hyp2f1"
407+
408+
def evaluate(self, x, mapping, component, index_values):
409+
"""Evaluate."""
410+
try:
411+
import scipy.special # type: ignore
412+
except ImportError:
413+
raise ValueError(
414+
"You must have scipy installed to evaluate hypergeometric functions in python."
415+
)
416+
func = getattr(scipy.special, self._name)
417+
return func(*(arg.evaluate(x, mapping, component, index_values) for arg in self.ufl_operands))
418+
419+
def __str__(self):
420+
"""Format as a string."""
421+
return f"{self._name}({self.ufl_operands[0]}, {self.ufl_operands[1]}, {self.ufl_operands[2]}, {self.ufl_operands[3]})"
422+
423+
424+
391425
@ufl_type(is_abstract=True, is_scalar=True, num_ops=2)
392426
class BesselFunction(Operator):
393427
"""Base class for all bessel functions."""

ufl/operators.py

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -51,6 +51,7 @@
5151
Cosh,
5252
Erf,
5353
Exp,
54+
Hypergeometric2F1,
5455
Ln,
5556
Sin,
5657
Sinh,

0 commit comments

Comments
 (0)