File: forms.py

package info (click to toggle)
fenics-dolfinx 1%3A0.10.0.post4-1exp1
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 6,028 kB
  • sloc: cpp: 36,535; python: 25,391; makefile: 226; sh: 171; xml: 55
file content (692 lines) | stat: -rw-r--r-- 25,195 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
# Copyright (C) 2017-2024 Chris N. Richardson, Garth N. Wells,
# Michal Habera and Jørgen S. Dokken
#
# This file is part of DOLFINx (https://www.fenicsproject.org)
#
# SPDX-License-Identifier:    LGPL-3.0-or-later

from __future__ import annotations

import collections
import types
import typing
from collections.abc import Sequence
from dataclasses import dataclass

from mpi4py import MPI

import numpy as np
import numpy.typing as npt

import ffcx
import ufl
from dolfinx import cpp as _cpp
from dolfinx import default_scalar_type, jit
from dolfinx.fem import IntegralType
from dolfinx.fem.function import Constant, Function, FunctionSpace

if typing.TYPE_CHECKING:
    # import dolfinx.mesh just when doing type checking to avoid
    # circular import
    from dolfinx.mesh import EntityMap as _EntityMap
    from dolfinx.mesh import Mesh, MeshTags


class Form:
    _cpp_object: (
        _cpp.fem.Form_complex64
        | _cpp.fem.Form_complex128
        | _cpp.fem.Form_float32
        | _cpp.fem.Form_float64
    )
    _code: str | list[str] | None

    def __init__(
        self,
        form: _cpp.fem.Form_complex64
        | _cpp.fem.Form_complex128
        | _cpp.fem.Form_float32
        | _cpp.fem.Form_float64,
        ufcx_form=None,
        code: str | list[str] | None = None,
        module: types.ModuleType | list[types.ModuleType] | None = None,
    ):
        """A finite element form.

        Note:
            Forms should normally be constructed using :func:`form` and
            not using this class initialiser. This class is combined
            with different base classes that depend on the scalar type
            used in the Form.

        Args:
            form: Compiled form object.
            ufcx_form: UFCx form.
            code: Form C++ code.
            module: CFFI module.
        """
        self._code = code
        self._ufcx_form = ufcx_form
        self._cpp_object = form
        self._module = module

    @property
    def ufcx_form(self):
        """The compiled ufcx_form object."""
        return self._ufcx_form

    @property
    def code(self) -> str | list[str] | None:
        """C code strings."""
        return self._code

    @property
    def module(self) -> types.ModuleType | list[types.ModuleType] | None:
        """The CFFI module"""
        return self._module

    @property
    def rank(self) -> int:
        return self._cpp_object.rank

    @property
    def function_spaces(self) -> list[FunctionSpace]:
        """Function spaces on which this form is defined."""
        return self._cpp_object.function_spaces

    @property
    def dtype(self) -> np.dtype:
        """Scalar type of this form."""
        return np.dtype(self._cpp_object.dtype)

    @property
    def mesh(self) -> _cpp.mesh.Mesh_float32 | _cpp.mesh.Mesh_float64:
        """Mesh on which this form is defined."""
        return self._cpp_object.mesh

    @property
    def integral_types(self):
        """Integral types in the form."""
        return self._cpp_object.integral_types

    def num_integrals(self, integral_type: IntegralType, kernel_index: int) -> int:
        """Number of integrals of a given type for a specific cell type.

        Args:
            integral_type: The type of integral to count.
            kernel_index: In the case of mixed topology, we have a kernel
                per cell type. For single-cell type meshes, this is zero.
        """
        return self._cpp_object.num_integrals(integral_type, kernel_index)


def get_integration_domains(
    integral_type: IntegralType,
    subdomain: MeshTags | list[tuple[int, np.ndarray]] | None,
    subdomain_ids: list[int],
) -> list[tuple[int, np.ndarray]]:
    """Get integration domains from subdomain data.

    The subdomain data is a MeshTags object consisting of markers, or
    ``None``. If it is ``None``, we do not pack any integration
    entities. Integration domains are defined as a list of tuples, where
    each input ``subdomain_ids`` is mapped to an array of integration
    entities, where an integration entity for a cell integral is the
    list of cells. For an exterior facet integral each integration
    entity is a tuple ``(cell_index, local_facet_index)``. For an
    interior facet integral each integration entity is a tuple
    ``(cell_index0, local_facet_index0, cell_index1,
    local_facet_index1)``. Where the first cell-facet pair is the
    ``'+'`` restriction, the second the ``'-'`` restriction.

    Args:
        integral_type: The type of integral to pack integration
            entities for.
        subdomain: A MeshTag with markers or manually specified
            integration domains.
        subdomain_ids: List of ids to integrate over.

    Returns:
        A list of entities to integrate over. For cell integrals, this is a
        list of cells. For exterior facet integrals, this is a list of
        (cell, local_facet) pairs. For interior facet integrals, this is a
        list of (cell0, local_facet0, cell1, local_facet1) tuples.
    """
    if subdomain is None:
        return []
    else:
        domains = []
        if not isinstance(subdomain, list):
            if integral_type in (IntegralType.exterior_facet, IntegralType.interior_facet):
                tdim = subdomain.topology.dim
                subdomain._cpp_object.topology.create_connectivity(tdim - 1, tdim)
                subdomain._cpp_object.topology.create_connectivity(tdim, tdim - 1)

            if integral_type is IntegralType.vertex:
                tdim = subdomain.topology.dim
                subdomain._cpp_object.topology.create_connectivity(0, tdim)
                subdomain._cpp_object.topology.create_connectivity(tdim, 0)

            if integral_type is IntegralType.ridge:
                tdim = subdomain.topology.dim
                subdomain._cpp_object.topology.create_connectivity(tdim - 2, tdim)
                subdomain._cpp_object.topology.create_connectivity(tdim, tdim - 2)

            # Compute integration domains only for each subdomain id in
            # the integrals. If a process has no integral entities,
            # insert an empty array.
            for id in subdomain_ids:
                integration_entities = _cpp.fem.compute_integration_domains(
                    integral_type,
                    subdomain._cpp_object.topology,
                    subdomain.find(id),
                )
                domains.append((id, integration_entities))
            return [(s[0], np.array(s[1])) for s in domains]
        else:
            return [(s[0], np.array(s[1])) for s in sorted(subdomain)]


def form_cpp_class(
    dtype: npt.DTypeLike,
) -> (
    _cpp.fem.Form_float32
    | _cpp.fem.Form_float64
    | _cpp.fem.Form_complex64
    | _cpp.fem.Form_complex128
):
    """Return the wrapped C++ class of a variational form of a specific
    scalar type.

    Args:
        dtype: Scalar type of the required form class.

    Returns:
        Wrapped C++ form class of the requested type.

    Note:
        This function is for advanced usage, typically when writing
        custom kernels using Numba or C.
    """
    if np.issubdtype(dtype, np.float32):
        return _cpp.fem.Form_float32
    elif np.issubdtype(dtype, np.float64):
        return _cpp.fem.Form_float64
    elif np.issubdtype(dtype, np.complex64):
        return _cpp.fem.Form_complex64
    elif np.issubdtype(dtype, np.complex128):
        return _cpp.fem.Form_complex128
    else:
        raise NotImplementedError(f"Type {dtype} not supported.")


_ufl_to_dolfinx_domain = {
    "cell": IntegralType.cell,
    "exterior_facet": IntegralType.exterior_facet,
    "interior_facet": IntegralType.interior_facet,
    "vertex": IntegralType.vertex,
    "ridge": IntegralType.ridge,
}


def mixed_topology_form(
    forms: Sequence[ufl.Form],
    dtype: npt.DTypeLike = default_scalar_type,
    form_compiler_options: dict | None = None,
    jit_options: dict | None = None,
    entity_maps: Sequence[_EntityMap] | None = None,
):
    """
    Create a mixed-topology from from an array of Forms.

    # FIXME: This function is a temporary hack for mixed-topology
    meshes. # It is needed because UFL does not know about
    mixed-topology meshes, # so we need to pass a list of forms for each
    cell type.

    Args:
        form: A list of UFL forms. Each form should be the same, just
            defined on different cell types.
        dtype: Scalar type to use for the compiled form.
        form_compiler_options: See :func:`ffcx_jit <dolfinx.jit.ffcx_jit>`
        jit_options: See :func:`ffcx_jit <dolfinx.jit.ffcx_jit>`.
        entity_maps: If any trial functions, test functions, or
            coefficients in the form are not defined over the same mesh
            as the integration domain (the domain associated with the
            measure), `entity_maps` must be supplied. For each mesh in
            the form, there should be an entity map relating entities in
            that mesh to the integration domain mesh.

    Returns:
        Compiled finite element Form.
    """

    if form_compiler_options is None:
        form_compiler_options = dict()

    form_compiler_options["scalar_type"] = dtype
    ftype = form_cpp_class(dtype)

    # Extract subdomain data from UFL form
    sd = next(iter(forms)).subdomain_data()
    (domain,) = list(sd.keys())  # Assuming single domain

    # Check that subdomain data for each integral type is the same
    for data in sd.get(domain).values():
        assert all([d is data[0] for d in data if d is not None])

    mesh = domain.ufl_cargo()

    ufcx_forms = []
    modules = []
    codes = []
    for form in forms:
        ufcx_form, module, code = jit.ffcx_jit(
            mesh.comm,
            form,
            form_compiler_options=form_compiler_options,
            jit_options=jit_options,
        )
        ufcx_forms.append(ufcx_form)
        modules.append(module)
        codes.append(code)

    # In a mixed-topology mesh, each form has the same C++ function
    # space, so we can extract it from any of them
    V = [arg.ufl_function_space()._cpp_object for arg in form.arguments()]

    # TODO coeffs, constants, subdomains, entity_maps
    f = ftype(
        [module.ffi.cast("uintptr_t", module.ffi.addressof(ufcx_form)) for ufcx_form in ufcx_forms],
        V,
        [],
        [],
        {},
        [],
        mesh,
    )
    return Form(f, ufcx_forms, codes, modules)


def form(
    form: ufl.Form | Sequence[ufl.Form] | Sequence[Sequence[ufl.Form]],
    dtype: npt.DTypeLike = default_scalar_type,
    form_compiler_options: dict | None = None,
    jit_options: dict | None = None,
    entity_maps: Sequence[_EntityMap] | None = None,
):
    """Create a Form or list of Forms.

    Args:
        form: A UFL form or iterable of UFL forms.
        dtype: Scalar type to use for the compiled form.
        form_compiler_options: See :func:`ffcx_jit <dolfinx.jit.ffcx_jit>`
        jit_options: See :func:`ffcx_jit <dolfinx.jit.ffcx_jit>`.
        entity_maps: If any trial functions, test functions, or
            coefficients in the form are not defined over the same mesh
            as the integration domain (the domain associated with the
            measure), `entity_maps` must be supplied. For each mesh in
            the form, there should be an entity map relating entities in
            that mesh to the integration domain mesh.

    Returns:
        Compiled finite element Form.

    Note:
        This function is responsible for the compilation of a UFL form
        (using FFCx) and attaching coefficients and domains specific
        data to the underlying C++ form. It dynamically create a
        :class:`Form` instance with an appropriate base class for the
        scalar type, e.g. :func:`_cpp.fem.Form_float64`.
    """
    if form_compiler_options is None:
        form_compiler_options = dict()

    form_compiler_options["scalar_type"] = dtype
    ftype = form_cpp_class(dtype)

    def _form(form):
        """Compile a single UFL form."""
        # Extract subdomain data from UFL form
        sd = form.subdomain_data()
        (domain,) = list(sd.keys())  # Assuming single domain

        # Check that subdomain data for each integral type is the same
        for data in sd.get(domain).values():
            assert all([d is data[0] for d in data if d is not None])

        msh = domain.ufl_cargo()
        if msh is None:
            raise RuntimeError("Expecting to find a Mesh in the form.")
        ufcx_form, module, code = jit.ffcx_jit(
            msh.comm, form, form_compiler_options=form_compiler_options, jit_options=jit_options
        )

        # For each argument in form extract its function space
        V = [arg.ufl_function_space()._cpp_object for arg in form.arguments()]
        part = form_compiler_options.get("part", "full")
        if part == "diagonal":
            V = [V[0]]

        # Prepare coefficients data. For every coefficient in form take
        # its C++ object.
        original_coeffs = form.coefficients()
        coeffs = [
            original_coeffs[ufcx_form.original_coefficient_positions[i]]._cpp_object
            for i in range(ufcx_form.num_coefficients)
        ]
        constants = [c._cpp_object for c in form.constants()]

        # Extract subdomain ids from ufcx_form
        subdomain_ids = {type: [] for type in sd.get(domain).keys()}
        integral_offsets = [ufcx_form.form_integral_offsets[i] for i in range(6)]
        for i in range(len(integral_offsets) - 1):
            integral_type = IntegralType(i)
            for j in range(integral_offsets[i], integral_offsets[i + 1]):
                subdomain_ids[integral_type.name].append(ufcx_form.form_integral_ids[j])

        # Subdomain markers (possibly empty list for some integral types)
        subdomains = {
            _ufl_to_dolfinx_domain[key]: get_integration_domains(
                _ufl_to_dolfinx_domain[key], subdomain_data[0], subdomain_ids[key]
            )
            for (key, subdomain_data) in sd.get(domain).items()
        }

        if entity_maps is None:
            _entity_maps = []
        else:
            _entity_maps = [entity_map._cpp_object for entity_map in entity_maps]

        f = ftype(
            [module.ffi.cast("uintptr_t", module.ffi.addressof(ufcx_form))],
            V,
            coeffs,
            constants,
            subdomains,
            _entity_maps,
            msh,
        )
        return Form(f, ufcx_form, code, module)

    def _zero_form(form):
        """Compile a single 'zero' UFL form, i.e. a form with no
        integrals."""
        V = [arg.ufl_function_space()._cpp_object for arg in form.arguments()]
        assert len(V) > 0
        msh = V[0].mesh

        f = ftype(
            spaces=V,
            integrals={},
            coefficients=[],
            constants=[],
            need_permutation_data=False,
            entity_maps=[],
            mesh=msh,
        )
        return Form(f)

    def _create_form(form):
        """Recursively convert ufl.Forms to dolfinx.fem.Form.

        Args:
            form: UFL form or list of UFL forms to extract DOLFINx forms
                from.

        Returns:
            A ``dolfinx.fem.Form`` or a list of ``dolfinx.fem.Form``.
        """
        if isinstance(form, ufl.Form):
            return _form(form)
        elif isinstance(form, ufl.ZeroBaseForm):
            return _zero_form(form)
        elif isinstance(form, collections.abc.Iterable):
            return list(map(lambda sub_form: _create_form(sub_form), form))
        else:
            return form

    return _create_form(form)


def extract_function_spaces(
    forms: Form | Sequence[Form] | Sequence[Sequence[Form]],
    index: int = 0,
) -> FunctionSpace | list[None | FunctionSpace]:
    """Extract common function spaces from an array of forms.

    If ``forms`` is a list of linear forms, this function returns of list
    of the corresponding test function spaces. If ``forms`` is a 2D
    array of bilinear forms, for ``index=0`` the list of common test
    function spaces for each row is returned, and if ``index=1`` the
    common trial function spaces for each column are returned.

    Args:
        forms: A list of forms or a 2D array of forms.
        index: Index of the function space to extract. If ``index=0``,
            the test function spaces are extracted, if ``index=1`` the
            trial function spaces are extracted.

    Returns:
        List of function spaces.
    """
    _forms = np.array(forms)
    if _forms.ndim == 0:
        form: Form = _forms.tolist()
        return form.function_spaces[0] if form is not None else None
    elif _forms.ndim == 1:
        assert index == 0, "Expected index=0 for 1D array of forms"
        for form in _forms:
            if form is not None:
                assert form.rank == 1, "Expected linear form"
        return [form.function_spaces[0] if form is not None else None for form in forms]  # type: ignore[union-attr]
    elif _forms.ndim == 2:
        assert index == 0 or index == 1, "Expected index=0 or index=1 for 2D array of forms"
        extract_spaces = np.vectorize(
            lambda form: form.function_spaces[index] if form is not None else None
        )
        V = extract_spaces(_forms)

        def unique_spaces(V):
            # Pick spaces from first column
            V0 = V[:, 0]

            # Iterate over each column
            for col in range(1, V.shape[1]):
                # Iterate over entry in column, updating if current
                # space is None, or where both spaces are not None check
                # that they are the same
                for row in range(V.shape[0]):
                    if V0[row] is None and V[row, col] is not None:
                        V0[row] = V[row, col]
                    elif V0[row] is not None and V[row, col] is not None:
                        assert V0[row] is V[row, col], "Cannot extract unique function spaces"
            return V0

        if index == 0:
            return list(unique_spaces(V))
        elif index == 1:
            return list(unique_spaces(V.transpose()))

    raise RuntimeError("Unsupported array of forms")


@dataclass
class CompiledForm:
    """Compiled UFL form without associated DOLFINx data."""

    ufl_form: ufl.Form  # The original ufl form
    ufcx_form: typing.Any  # The compiled form
    module: typing.Any  #  The module
    code: str  # The source code
    dtype: npt.DTypeLike  # data type used for the `ufcx_form`


def compile_form(
    comm: MPI.Intracomm,
    form: ufl.Form,
    form_compiler_options: dict | None = {"scalar_type": default_scalar_type},
    jit_options: dict | None = None,
) -> CompiledForm:
    """Compile UFL form without associated DOLFINx data.

    Args:
        comm: The MPI communicator used when compiling the form
        form: The UFL form to compile
        form_compiler_options: See :func:`ffcx_jit <dolfinx.jit.ffcx_jit>`
        jit_options: See :func:`ffcx_jit <dolfinx.jit.ffcx_jit>`.
    """
    p_ffcx = ffcx.get_options(form_compiler_options)
    p_jit = jit.get_options(jit_options)
    ufcx_form, module, code = jit.ffcx_jit(comm, form, p_ffcx, p_jit)
    scalar_type: npt.DTypeLike = p_ffcx["scalar_type"]  # type: ignore [assignment]
    return CompiledForm(form, ufcx_form, module, code, scalar_type)


def form_cpp_creator(
    dtype: npt.DTypeLike,
) -> (
    _cpp.fem.Form_float32
    | _cpp.fem.Form_float64
    | _cpp.fem.Form_complex64
    | _cpp.fem.Form_complex128
):
    """Return the wrapped C++ constructor for creating a variational form
    of a specific scalar type.

    Args:
        dtype: Scalar type of the required form class.

    Returns:
        Wrapped C++ form class of the requested type.

    Note:
        This function is for advanced usage, typically when writing
        custom kernels using Numba or C.
    """
    if np.issubdtype(dtype, np.float32):
        return _cpp.fem.create_form_float32
    elif np.issubdtype(dtype, np.float64):
        return _cpp.fem.create_form_float64
    elif np.issubdtype(dtype, np.complex64):
        return _cpp.fem.create_form_complex64
    elif np.issubdtype(dtype, np.complex128):
        return _cpp.fem.create_form_complex128
    else:
        raise NotImplementedError(f"Type {dtype} not supported.")


def create_form(
    form: CompiledForm,
    function_spaces: list[FunctionSpace],
    msh: Mesh,
    subdomains: dict[IntegralType, list[tuple[int, np.ndarray]]],
    coefficient_map: dict[ufl.Coefficient, Function],
    constant_map: dict[ufl.Constant, Constant],
    entity_maps: Sequence[_EntityMap] | None = None,
) -> Form:
    """Create a Form object from a data-independent compiled form.

    Args:
        form: Compiled ufl form function_spaces: List of function spaces
        associated with the
            form. Should match the number of arguments in the form.
        msh: Mesh to associate form with. subdomains: A map from
        integral type to a list of pairs, where
            each pair corresponds to a subdomain id and the set of of
            integration entities to integrate over. Can be computed with
            {py:func}`dolfinx.fem.compute_integration_domains`.
        coefficient_map: Map from UFL coefficient to function with data.
        constant_map: Map from UFL constant to constant with data.
            to the integration domain ``msh``. The value of the map is
            an array of integers, where the i-th entry is the entity in
            the key mesh.
        entity_maps: Entity maps to support cases where forms involve
            sub-meshes.

    Return:
        A Form object.
    """
    if entity_maps is None:
        _entity_maps = []
    else:
        _entity_maps = [entity_map._cpp_object for entity_map in entity_maps]

    _subdomain_data = subdomains.copy()
    for _, idomain in _subdomain_data.items():
        idomain.sort(key=lambda x: x[0])

    # Extract all coefficients of the compiled form in correct order
    coefficients = {}
    original_coefficients = ufl.algorithms.extract_coefficients(form.ufl_form)
    num_coefficients = form.ufcx_form.num_coefficients
    for c in range(num_coefficients):
        original_index = form.ufcx_form.original_coefficient_positions[c]
        original_coeff = original_coefficients[original_index]
        try:
            coefficients[f"w{c}"] = coefficient_map[original_coeff]._cpp_object
        except KeyError:
            raise RuntimeError(f"Missing coefficient {original_coeff}")

    # Extract all constants of the compiled form in correct order
    # NOTE: Constants are not eliminated
    original_constants = ufl.algorithms.analysis.extract_constants(form.ufl_form)
    num_constants = form.ufcx_form.num_constants
    if num_constants != len(original_constants):
        raise RuntimeError(
            f"Number of constants in compiled form ({num_constants})",
            f"does not match the original form {len(original_constants)}",
        )
    constants = {}
    for counter, constant in enumerate(original_constants):
        try:
            mapped_constant = constant_map[constant]
            constants[f"c{counter}"] = mapped_constant._cpp_object
        except KeyError:
            raise RuntimeError(f"Missing constant {constant}")

    ftype = form_cpp_creator(form.dtype)
    f = ftype(
        form.module.ffi.cast("uintptr_t", form.module.ffi.addressof(form.ufcx_form)),
        [fs._cpp_object for fs in function_spaces],
        coefficients,
        constants,
        _subdomain_data,
        _entity_maps,
        msh._cpp_object,
    )
    return Form(f, form.ufcx_form, form.code)


def derivative_block(
    F: ufl.Form | Sequence[ufl.Form],
    u: Function | Sequence[Function],
    du: ufl.Argument | Sequence[ufl.Argument] | None = None,
) -> ufl.Form | Sequence[Sequence[ufl.Form]]:
    """Return the UFL derivative of a (list of) UFL rank one form(s).

    This is commonly used to derive a block Jacobian from a block
    residual.

    If `F_i` is a list of forms, the Jacobian is a list of lists with
    :math:`J_{ij} = \\frac{\\partial F_i}{u_j}[\\delta u_j]` using
    `ufl.derivative` called component-wise.

    If `F` is a form, the Jacobian is computed as :math:`J =
    \\frac{\\partial F}{\\partial u}[\\delta u]`. This is identical to
    calling `ufl.derivative` directly.
    """
    if isinstance(F, ufl.Form):
        if not isinstance(u, Function):
            raise ValueError("Must provide a single function when F is a UFL form")
        if du is None:
            du = ufl.TrialFunction(u.function_space)
        return ufl.derivative(F, u, du)
    else:
        assert all([isinstance(Fi, ufl.Form) for Fi in F]), "F must be a sequence of UFL forms"
        assert len(F) == len(u), "Number of forms and functions must be equal"
        if du is not None:
            assert len(F) == len(du), "Number of forms and du must be equal"
        else:
            du = [ufl.TrialFunction(u_i.function_space) for u_i in u]
        return [[ufl.derivative(Fi, u_j, du_j) for u_j, du_j in zip(u, du)] for Fi in F]