-
Notifications
You must be signed in to change notification settings - Fork 161
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
Apply Diriclet boundary on the Cofunction RHS. #3754
Draft
Ig-dolci
wants to merge
1
commit into
master
Choose a base branch
from
dolci/fix_issue3498
base: master
Could not load branches
Branch not found: {{ refName }}
Loading
Could not load tags
Nothing to show
Loading
Are you sure you want to change the base?
Some commits from the old base branch may be removed from the timeline,
and old review comments may become outdated.
Draft
Changes from all commits
Commits
File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -10,6 +10,7 @@ | |
DEFAULT_SNES_PARAMETERS | ||
) | ||
from firedrake.function import Function | ||
from firedrake.cofunction import Cofunction | ||
from firedrake.functionspace import RestrictedFunctionSpace | ||
from firedrake.ufl_expr import TrialFunction, TestFunction | ||
from firedrake.bcs import DirichletBC, EquationBC | ||
|
@@ -305,6 +306,10 @@ def solve(self, bounds=None): | |
|
||
for dbc in problem.dirichlet_bcs(): | ||
dbc.apply(problem.u_restrict) | ||
for coeff in coefficients: | ||
if isinstance(coeff, Cofunction): | ||
# Apply the DirichletBC to the right hand side of the equation. | ||
dbc.apply(coeff) | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. This will have side effects on
|
||
|
||
if bounds is not None: | ||
lower, upper = bounds | ||
|
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -60,3 +60,43 @@ def test_scalar_cofunction_zero_with_subset(V): | |
assert f is g | ||
assert np.allclose(f.dat.data_ro[:2], 0.0) | ||
assert np.allclose(f.dat.data_ro[2:], 1.0) | ||
|
||
|
||
def test_diriclet_bc_rhs(V): | ||
# Issue https://github.com/firedrakeproject/firedrake/issues/3498 | ||
# Apply DirichletBC to RHS (Cofunction) in LinearVariationalSolver | ||
mesh = UnitIntervalMesh(2) | ||
space = FunctionSpace(mesh, "Lagrange", 1) | ||
test, trial = TestFunction(space), TrialFunction(space) | ||
|
||
# Form RHS | ||
u = Function(space, name="u") | ||
problem = LinearVariationalProblem( | ||
inner(trial, test) * dx, inner(Constant(1.0), test) * dx, u, | ||
DirichletBC(space, 0.0, "on_boundary")) | ||
solver = LinearVariationalSolver(problem) | ||
solver.solve() | ||
|
||
assert np.allclose(assemble(inner(u, u) * ds), 0.0) | ||
|
||
# Cofunction RHS | ||
b = assemble(inner(Constant(1.0), test) * dx) | ||
u = Function(space, name="u") | ||
problem = LinearVariationalProblem( | ||
inner(trial, test) * dx, b, u, | ||
DirichletBC(space, 0.0, "on_boundary")) | ||
solver = LinearVariationalSolver(problem) | ||
solver.solve() | ||
|
||
assert np.allclose(assemble(inner(u, u) * ds), 0.0) | ||
|
||
# FormSum RHS | ||
b = assemble(inner(Constant(0.5), test) * dx) + inner(Constant(0.5), test) * dx | ||
u = Function(space, name="u") | ||
problem = LinearVariationalProblem( | ||
inner(trial, test) * dx, b, u, | ||
DirichletBC(space, 0.0, "on_boundary")) | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. The test should have non-homogeneous bcs |
||
solver = LinearVariationalSolver(problem) | ||
solver.solve() | ||
|
||
assert np.allclose(assemble(inner(u, u) * ds), 0.0) |
Oops, something went wrong.
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
This does not look right. The residual form might contain many cofunctions in function spaces different from the one of the solution. This should be done before calling the solver, and only on the Cofunction RHS. Also it is likely that you need dbc.zero() instead of dbc.apply(). Also see the comment I left in #3662