Skip to content

Commit 7e1d8ad

Browse files
authored
Merge pull request #4719 from firedrakeproject/connorjward/merge-release
Merge release into main
2 parents 497b0e0 + 25e7998 commit 7e1d8ad

File tree

19 files changed

+293
-199
lines changed

19 files changed

+293
-199
lines changed

AUTHORS.rst

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -64,6 +64,8 @@ Cyrus Cheng
6464

6565
Teodoro Fields Collin
6666

67+
Leo Collins
68+
6769
Colin J. Cotter...............<https://www.imperial.ac.uk/people/colin.cotter>
6870

6971
Joshua Coutinho

demos/saddle_point_pc/saddle_point_systems.py.rst

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -449,10 +449,10 @@ we then use to solve the problem, we can provide one ourselves.
449449
Recall that :math:`S` is spectrally a Laplacian only in a
450450
discontinuous space. A natural choice is therefore to use an interior
451451
penalty DG formulation for the Laplacian term on the block of the scalar
452-
variable. We can provide it as an :class:`~.AuxiliaryOperatorPC` via a python preconditioner. ::
452+
variable. We can provide it as an :class:`~.AuxiliaryOperatorPC` via a python preconditioner. Note that the ```form``` method in ```AuxiliaryOperatorPC``` takes the test functions as the first argument and the trial functions as the second argument, which is the reverse of the usual convention. ::
453453

454454
class DGLaplacian(AuxiliaryOperatorPC):
455-
def form(self, pc, u, v):
455+
def form(self, pc, v, u):
456456
W = u.function_space()
457457
n = FacetNormal(W.mesh())
458458
alpha = Constant(4.0)

docs/source/contact.rst

Lines changed: 0 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -29,14 +29,3 @@ This is open to all, but you must `request an invite to join
2929
the channel. Should this not work for whatever reason, please get in
3030
touch via `GitHub discussions
3131
<https://github.com/firedrakeproject/firedrake/discussions>`__.
32-
33-
34-
Mailing list
35-
------------
36-
37-
Please join the Firedrake mailing list: [email protected]. This is a
38-
very low traffic list but it does carry important announcements, for example
39-
when we change a user-facing interface. Join the list on `this page
40-
<mailing_list_>`_.
41-
42-
.. _mailing_list: https://mailman.ic.ac.uk/mailman/listinfo/firedrake

docs/source/images/leocollins.jpeg

23.6 KB
Loading

docs/source/install.rst

Lines changed: 14 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -426,6 +426,20 @@ To install Firedrake with SLEPc support you should:
426426

427427
$ pip install --no-binary h5py 'firedrake[check,slepc]'
428428

429+
VTK
430+
~~~
431+
432+
To install Firedrake with VTK, it should be installed using the ``vtk`` optional
433+
dependency. For example::
434+
435+
$ pip install --no-binary h5py 'firedrake[check,vtk]'
436+
437+
.. warning::
438+
439+
VTK make releases sporadically so will not always support the latest version
440+
of Python. This is commonly an issue on macOS where homebrew will use the
441+
latest Python very soon after it is released. To fix this you should use
442+
an older version of Python (e.g. 3.13 instead of 3.14).
429443

430444
PyTorch
431445
~~~~~~~

docs/source/team.ini

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -49,6 +49,7 @@ India Marsden: https://www.maths.ox.ac.uk/people/india.marsden
4949
Daiane I. Dolci: https://www.imperial.ac.uk/people/d.dolci
5050
Joshua Hope-Collins: https://www.imperial.ac.uk/people/joshua.hope-collins13
5151
Umberto Zerbinati: https://www.uzerbinati.eu
52+
Leo Collins:
5253

5354
[inactive-team]
5455
Lawrence Mitchell: https://www.wence.uk/

firedrake/adjoint_utils/function.py

Lines changed: 6 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -278,7 +278,12 @@ def _ad_dot(self, other, options=None):
278278
def _ad_assign_numpy(dst, src, offset):
279279
range_begin, range_end = dst.dat.dataset.layout_vec.getOwnershipRange()
280280
m_a_local = src[offset + range_begin:offset + range_end]
281-
dst.dat.data_wo[...] = m_a_local.reshape(dst.dat.data_wo.shape)
281+
if dst.function_space().ufl_element().family() == "Real":
282+
# Real space keeps a redundant copy of the data on every rank
283+
comm = dst.function_space().mesh()._comm
284+
dst.dat.data_wo[...] = comm.bcast(m_a_local, root=0)
285+
else:
286+
dst.dat.data_wo[...] = m_a_local.reshape(dst.dat.data_wo.shape)
282287
offset += dst.dat.dataset.layout_vec.size
283288
return dst, offset
284289

firedrake/function.py

Lines changed: 6 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -563,8 +563,12 @@ def evaluate(self, coord, mapping, component, index_values):
563563
# Called by UFL when evaluating expressions at coordinates
564564
if component or index_values:
565565
raise NotImplementedError("Unsupported arguments when attempting to evaluate Function.")
566+
coord = np.asarray(coord, dtype=utils.ScalarType)
566567
evaluator = PointEvaluator(self.function_space().mesh(), coord)
567-
return evaluator.evaluate(self)
568+
result = evaluator.evaluate(self)
569+
if len(coord.shape) == 1:
570+
result = result.squeeze(axis=0)
571+
return result
568572

569573
def at(self, arg, *args, **kwargs):
570574
warnings.warn(
@@ -827,7 +831,7 @@ def evaluate(self, function: Function) -> np.ndarray | Tuple[np.ndarray, ...]:
827831
f_at_points = assemble(interpolate(function, P0DG))
828832
f_at_points_io = Function(P0DG_io).assign(np.nan)
829833
f_at_points_io.interpolate(f_at_points)
830-
result = f_at_points_io.dat.data_ro
834+
result = f_at_points_io.dat.data_ro.copy()
831835

832836
# If redundant, all points are now on rank 0, so we broadcast the result
833837
if self.redundant and self.mesh.comm.size > 1:

firedrake/output/vtk_output.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -347,7 +347,7 @@ def get_array(function):
347347
return array
348348

349349

350-
class VTKFile(object):
350+
class VTKFile:
351351
_header = (b'<?xml version="1.0" ?>\n'
352352
b'<VTKFile type="Collection" version="0.1" '
353353
b'byte_order="LittleEndian">\n'

firedrake/preconditioners/hypre_ads.py

Lines changed: 24 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -1,11 +1,12 @@
11
from firedrake.preconditioners.base import PCBase
2+
from firedrake.preconditioners.fdm import tabulate_exterior_derivative
23
from firedrake.petsc import PETSc
34
from firedrake.function import Function
45
from firedrake.ufl_expr import TrialFunction
56
from firedrake.dmhooks import get_function_space
67
from firedrake.preconditioners.hypre_ams import chop
78
from firedrake.interpolation import interpolate
8-
from finat.ufl import VectorElement
9+
from finat.ufl import FiniteElement, TensorElement, VectorElement
910
from ufl import grad, curl, SpatialCoordinate
1011
from pyop2.utils import as_tuple
1112

@@ -27,16 +28,34 @@ def initialize(self, obj):
2728
if formdegree != 2 or degree != 1:
2829
raise ValueError("Hypre ADS requires lowest order RT elements! (not %s of degree %d)" % (family, degree))
2930

30-
P1 = V.reconstruct(family="Lagrange", degree=1)
31-
NC1 = V.reconstruct(family="N1curl" if mesh.ufl_cell().is_simplex else "NCE", degree=1)
31+
# Get the auxiliary Nedelec and Lagrange spaces and the coordinate space
32+
cell = V.ufl_element().cell
33+
NC1_element = FiniteElement("N1curl" if cell.is_simplex else "NCE", cell=cell, degree=1)
34+
P1_element = FiniteElement("Lagrange", cell=cell, degree=1)
35+
coords_element = VectorElement(P1_element, dim=mesh.geometric_dimension)
36+
if V.shape:
37+
NC1_element = TensorElement(NC1_element, shape=V.shape)
38+
P1_element = TensorElement(P1_element, shape=V.shape)
39+
40+
NC1 = V.reconstruct(element=NC1_element)
41+
P1 = V.reconstruct(element=P1_element)
42+
VectorP1 = V.reconstruct(element=coords_element)
43+
3244
G_callback = appctx.get("get_gradient", None)
3345
if G_callback is None:
34-
G = chop(assemble(interpolate(grad(TrialFunction(P1)), NC1)).petscmat)
46+
try:
47+
G = chop(assemble(interpolate(grad(TrialFunction(P1)), NC1)).petscmat)
48+
except NotImplementedError:
49+
# dual evaluation not yet implemented see https://github.com/firedrakeproject/fiat/issues/109
50+
G = tabulate_exterior_derivative(P1, NC1)
3551
else:
3652
G = G_callback(P1, NC1)
3753
C_callback = appctx.get("get_curl", None)
3854
if C_callback is None:
39-
C = chop(assemble(interpolate(curl(TrialFunction(NC1)), V)).petscmat)
55+
try:
56+
C = chop(assemble(interpolate(curl(TrialFunction(NC1)), V)).petscmat)
57+
except NotImplementedError:
58+
C = tabulate_exterior_derivative(NC1, V)
4059
else:
4160
C = C_callback(NC1, V)
4261

@@ -50,7 +69,6 @@ def initialize(self, obj):
5069
pc.setHYPREDiscreteGradient(G)
5170
pc.setHYPREDiscreteCurl(C)
5271

53-
VectorP1 = P1.reconstruct(element=VectorElement(P1.ufl_element()))
5472
coords = Function(VectorP1).interpolate(SpatialCoordinate(mesh))
5573
pc.setCoordinates(coords.dat.data_ro.copy())
5674

0 commit comments

Comments
 (0)