-
-
Notifications
You must be signed in to change notification settings - Fork 75
Open
Description
The following MWE, based on the post:
https://fenicsproject.discourse.group/t/shape-optimization-in-complex-dolfinx/18345/5
yields a zero-derivative for complex mode.
import ufl
from utils import LagrangeElement
import numpy as np
scalar_type=np.float64
mesh = ufl.Mesh(LagrangeElement(ufl.triangle, 1, shape=(2, )))
X = ufl.SpatialCoordinate(mesh)
W = ufl.FunctionSpace(mesh, mesh.ufl_coordinate_element())
volume = ufl.Constant(mesh)*ufl.dx
ddvol = ufl.derivative(volume, X, ufl.conj(ufl.TestFunction(W)))
complex_mode = np.issubdtype(scalar_type, np.complexfloating)
fd = ufl.algorithms.compute_form_data(ddvol,
do_apply_function_pullbacks=True,
do_apply_integral_scaling=True,
do_apply_geometry_lowering=True,
preserve_geometry_types=(ufl.classes.Jacobian,),
do_apply_restrictions=True,
do_append_everywhere_integrals=False,
complex_mode=complex_mode)
print(complex_mode, str(fd.preprocessed_form))
scalar_type=np.complex64
complex_mode = np.issubdtype(scalar_type, np.complexfloating)
fd = ufl.algorithms.compute_form_data(ddvol,
do_apply_function_pullbacks=True,
do_apply_integral_scaling=True,
do_apply_geometry_lowering=True,
preserve_geometry_types=(ufl.classes.Jacobian,),
do_apply_restrictions=True,
do_append_everywhere_integrals=False,
complex_mode=complex_mode)
print(complex_mode, str(fd.preprocessed_form))output:
False { c_0 * weight * ((reference_grad(reference_value(v_0)))[0, 0] * J[1, 1] + J[0, 0] * (reference_grad(reference_value(v_0)))[1, 1] + -1 * ((reference_grad(reference_value(v_0)))[0, 1] * J[1, 0] + J[0, 1] * (reference_grad(reference_value(v_0)))[1, 0])) * (((J[0, 0] * J[1, 1] + -1 * J[0, 1] * J[1, 0]) == (0)) ? (0) : (((J[0, 0] * J[1, 1] + -1 * J[0, 1] * J[1, 0]) < (0)) ? (-1) : (1))) } * dx(<Mesh #0>[('otherwise',)], {'estimated_polynomial_degree': 1})
True <empty Form>Metadata
Metadata
Assignees
Labels
No labels