Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add generic nodal/modal bilinear form routines #103

Open
wants to merge 19 commits into
base: main
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from 7 commits
Commits
Show all changes
19 commits
Select commit Hold shift + click to select a range
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
8 changes: 4 additions & 4 deletions contrib/weights-from-festa-sommariva.py
Original file line number Diff line number Diff line change
Expand Up @@ -76,17 +76,17 @@ def generate_festa_sommariva_quadrature_rules(outfile):
#
# DEGREE: <degree>

_re_degree = re.compile(r"DEGREE:\s+(\d+)")
_re_xw = re.compile(r"xw\s?=\s?\[(.+?)\];", re.DOTALL)
re_degree = re.compile(r"DEGREE:\s+(\d+)")
re_xw = re.compile(r"xw\s?=\s?\[(.+?)\];", re.DOTALL)

# NOTE: some degrees have multiple quadrature rules with a repeated
# header, so we just take the unique ones here
degrees = np.unique(
np.fromiter(_re_degree.findall(mfile), dtype=np.int64)
np.fromiter(re_degree.findall(mfile), dtype=np.int64)
)

rules = {}
for imatch, match in enumerate(_re_xw.findall(mfile)):
for imatch, match in enumerate(re_xw.findall(mfile)):
d = degrees[imatch]
assert d == imatch + 1, (d, imatch + 1)

Expand Down
4 changes: 4 additions & 0 deletions modepy/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -32,10 +32,12 @@
inverse_mass_matrix,
mass_matrix,
modal_mass_matrix_for_face,
modal_quad_bilinear_form,
modal_quad_mass_matrix,
modal_quad_mass_matrix_for_face,
multi_vandermonde,
nodal_mass_matrix_for_face,
nodal_quad_bilinear_form,
nodal_quad_mass_matrix,
nodal_quad_mass_matrix_for_face,
resampling_matrix,
Expand Down Expand Up @@ -158,11 +160,13 @@
"legendre_gauss_tensor_product_nodes",
"mass_matrix",
"modal_mass_matrix_for_face",
"modal_quad_bilinear_form",
"modal_quad_mass_matrix",
"modal_quad_mass_matrix_for_face",
"monomial_basis_for_space",
"multi_vandermonde",
"nodal_mass_matrix_for_face",
"nodal_quad_bilinear_form",
"nodal_quad_mass_matrix",
"nodal_quad_mass_matrix_for_face",
"node_tuples_for_space",
Expand Down
85 changes: 71 additions & 14 deletions modepy/matrices.py
Copy link
Owner

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

While we're at it, deprecate the *face_mass* things, too, and (maybe) have a basis-with-node-transform-composer?

Original file line number Diff line number Diff line change
Expand Up @@ -66,6 +66,8 @@
.. autofunction:: modal_quad_mass_matrix_for_face
.. autofunction:: nodal_quad_mass_matrix_for_face
.. autofunction:: spectral_diag_nodal_mass_matrix
.. autofunction:: modal_quad_bilinear_form
.. autofunction:: nodal_quad_bilinear_form

Differentiation is also convenient to express by using :math:`V^{-1}` to
obtain modal values and then using a Vandermonde matrix for the derivatives
Expand Down Expand Up @@ -351,6 +353,68 @@ def mass_matrix(
return la.inv(inverse_mass_matrix(basis, nodes))


def modal_quad_bilinear_form(
Copy link
Owner

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I'm thinking about the naming of this. I'm not sure "bilinear form" is that great a fit, because the trial functions don't appear. I thought of "projection", but that also misses the mark, because the test functions could be derivatives. In a word, I'm not sure. I can probably be convinced to go along with "bilinear form", but I figured it might be good to at least discuss the naming.

(@alexfikl, if you're reading along, please feel free to chime in too)

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I'm also not the biggest fan of the name, but I couldn't come up with anything better. modal_quadrature_operator maybe?

Copy link
Collaborator

@alexfikl alexfikl Dec 14, 2024

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Where is this modal_quad_bilinear_form supposed to be used? I'm guessing the nodal_quad_bilinear_form is used to define a proper bilinear form somewhere?

EDIT: I mostly meant "proper" in the sense of defining an inner product, but maybe that's not super relevant..?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

modal_quad_bilinear_form doesn't have any uses right now other than returning the operator that's transformed into a nodal operator.

nodal_quad_bilinear_form will be used in a grudge update that evaluates operators using quadrature (inducer/grudge#362 is the staging area for this, although a little outdated at this point). While nodal_quadrature_operator could be used for everything, nodal_quad_bilinear_form makes the corresponding grudge code cleaner when there is no overintegration.

TLDR; modal_quad_bilinear_form's current use is in nodal_quad_bilinear_form. nodal_quad_bilinear_form is used to evaluate DG operators in grudge when there is no overintegration. nodal_quadrature_operator is used when there is overintegration.

Also, the original comment is outdated now since modal_quad_bilinear_form does use trial functions (as of c063d23).

Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Thanks for explaining! I think I got a better idea of what this is for.. but no better recommendation for a name :\

The only thing that comes to mind is that if modal_quad_bilinear_form is not actually supposed to be used outside of modepy at the moment, we can just inline it in nodal_quad_bilinear_form?

quadrature: Quadrature,
test_functions: Sequence[Callable[[np.ndarray], np.ndarray]]
) -> np.ndarray:
r"""Using the *quadrature*, provide a matrix :math:`A` that
satisfies:

.. math::

\displaystyle (A \boldsymbol u)_i = \sum_j w_j \phi_i(r_j) u_j,

where :math:`\phi_i` are the *test_functions* at the nodes
:math:`r_j` of *quadrature*, with corresponding weights :math:`w_j`.
"""
modal_operator = np.empty((len(test_functions), len(quadrature.weights)))

for i, test_f in enumerate(test_functions):
modal_operator[i] = test_f(quadrature.nodes) * quadrature.weights

return modal_operator
Comment on lines +373 to +378
Copy link
Owner

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

return np.array([...])?



def nodal_quad_bilinear_form(
test_functions: Sequence[Callable[[np.ndarray], np.ndarray]],
trial_functions: Sequence[Callable[[np.ndarray], np.ndarray]],
quadrature: Quadrature,
nodes: np.ndarray,
uses_quadrature_domain: bool = False
) -> np.ndarray:
a-alveyblanc marked this conversation as resolved.
Show resolved Hide resolved
r"""Using *quadrature*, provide a matrix :math:`A` that satisfies:

.. math::

\displaystyle (A \boldsymbol u)_i = \sum_j w_j \phi_i(r_j) u_j,

where :math:`\phi_i` are the Lagrange basis functions obtained from
*test_functions* at *nodes*, :math:`w_j` and :math:`r_j` are the weights
and nodes from *quadrature*, and :math:`u_j` are point values of the trial
function at either *nodes* or the *quadrature* nodes depending on whether
a quadrature domain is used (as signified by *uses_quadrature_domain*).

If *uses_quadrature_domain* is set to False, then an interpolation operator
is used to make the operator :math:`N \times N` where :math:`N` is the
number of nodes in *nodes*. Otherwise, the operator has shape :math:`N
\times N_q` where :math:`N_q` is the number of nodes associated with
*quadrature*. Default behavior is to assume a quadrature domain is not used.
"""
if len(test_functions) != nodes.shape[1]:
raise ValueError("volume_nodes not unisolvent with trial_functions")

vdm_out = vandermonde(trial_functions, nodes)
a-alveyblanc marked this conversation as resolved.
Show resolved Hide resolved

nodal_operator = la.solve(
vdm_out.T, modal_quad_bilinear_form(quadrature, test_functions))

if uses_quadrature_domain:
return nodal_operator
else:
return nodal_operator @ resampling_matrix(
trial_functions, quadrature.nodes, nodes)


def modal_quad_mass_matrix(
Copy link
Owner

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I'd vote in favor of deprecating the _mass_matrix functions.

quadrature: Quadrature,
test_functions: Sequence[Callable[[np.ndarray], np.ndarray]],
Expand All @@ -367,12 +431,7 @@ def modal_quad_mass_matrix(

.. versionadded :: 2024.2
"""
modal_mass_matrix = np.empty((len(test_functions), len(quadrature.weights)))

for i, test_f in enumerate(test_functions):
modal_mass_matrix[i] = test_f(quadrature.nodes) * quadrature.weights

return modal_mass_matrix
return modal_quad_bilinear_form(quadrature, test_functions)


def nodal_quad_mass_matrix(
Expand All @@ -399,13 +458,12 @@ def nodal_quad_mass_matrix(
if nodes is None:
nodes = quadrature.nodes

if len(test_functions) != nodes.shape[1]:
raise ValueError("volume_nodes not unisolvent with test_functions")

vdm = vandermonde(test_functions, nodes)

return la.solve(vdm.T,
modal_quad_mass_matrix(quadrature, test_functions))
return nodal_quad_bilinear_form(
test_functions,
test_functions,
quadrature,
nodes,
uses_quadrature_domain=True)


def spectral_diag_nodal_mass_matrix(
Expand Down Expand Up @@ -549,5 +607,4 @@ def nodal_quad_mass_matrix_for_face(
return la.solve(vol_vdm.T,
modal_quad_mass_matrix_for_face(face, face_quad, test_functions))


# vim: foldmethod=marker
64 changes: 64 additions & 0 deletions modepy/test/test_matrices.py
Original file line number Diff line number Diff line change
Expand Up @@ -118,6 +118,70 @@ def test_tensor_product_diag_mass_matrix(shape: mp.Shape) -> None:
assert np.max(err) < 1e-14


@pytest.mark.parametrize("shape_cls", [mp.Hypercube, mp.Simplex])
@pytest.mark.parametrize("dim", [1, 2, 3])
@pytest.mark.parametrize("order", [0, 1, 2, 4])
@pytest.mark.parametrize("nodes_on_bdry", [False, True])
@pytest.mark.parametrize("test_weak_d_dr", [False, True])
def test_bilinear_forms(
shape_cls: type[mp.Shape],
dim: int,
order: int,
nodes_on_bdry: bool,
test_weak_d_dr: bool
) -> None:
shape = shape_cls(dim)
space = mp.space_for_shape(shape, order)

quad_space = mp.space_for_shape(shape, 2*order)
quad = mp.quadrature_for_space(quad_space, shape)

if isinstance(shape, mp.Hypercube) or shape == mp.Simplex(1):
if nodes_on_bdry:
nodes = mp.legendre_gauss_lobatto_tensor_product_nodes(
shape.dim, order,
)
else:
nodes = mp.legendre_gauss_tensor_product_nodes(shape.dim, order)
elif isinstance(shape, mp.Simplex):
if nodes_on_bdry:
nodes = mp.warp_and_blend_nodes(shape.dim, order)
else:
nodes = mp.VioreanuRokhlinSimplexQuadrature(order, shape.dim).nodes
else:
raise AssertionError()
a-alveyblanc marked this conversation as resolved.
Show resolved Hide resolved

basis = mp.orthonormal_basis_for_space(space, shape)

if test_weak_d_dr and order not in [0, 1]:
mass_inv = mp.inverse_mass_matrix(basis, nodes)

for ax in range(dim):
f = 1 - nodes[ax]**2
fp = -2*nodes[ax]

weak_operator = mp.nodal_quad_bilinear_form(
basis.derivatives(ax),
basis.functions,
quad,
nodes,
uses_quadrature_domain=False)

err = la.norm(mass_inv @ weak_operator.T @ f - fp) / la.norm(fp)
assert err <= 1e-12
else:
quad_mass_mat = mp.nodal_quad_bilinear_form(
basis.functions,
basis.functions,
quad,
nodes,
uses_quadrature_domain=False)

vdm_mass_mat = mp.mass_matrix(basis, nodes)
err = la.norm(quad_mass_mat - vdm_mass_mat) / la.norm(vdm_mass_mat)
assert err < 1e-14


# You can test individual routines by typing
# $ python test_modes.py 'test_routine()'

Expand Down
Loading