-
Notifications
You must be signed in to change notification settings - Fork 126
Use Jacobian-free Newton-Krylov for SDIRK time integrators #2505
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
Open
DanielDoehring
wants to merge
35
commits into
trixi-framework:main
Choose a base branch
from
DanielDoehring:ExamplesNewtonKrylov
base: main
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.
Open
Changes from 19 commits
Commits
Show all changes
35 commits
Select commit
Hold shift + click to select a range
c4cbb5a
example linear newton krylov
DanielDoehring 8a4fdcd
Nonlinear example
DanielDoehring 94f250a
comment
DanielDoehring 2a51e91
Update test/Project.toml
DanielDoehring 16dc218
Update test/Project.toml
DanielDoehring 238afe3
check different sciml
DanielDoehring 1221c94
v
DanielDoehring 3453202
v
DanielDoehring 7512dce
v
DanielDoehring 57e5989
v
DanielDoehring 1ae7504
v
DanielDoehring a1ffd0d
v
DanielDoehring 0adb284
v
DanielDoehring b180056
v
DanielDoehring 83dce58
v
DanielDoehring 41f0b7c
v
DanielDoehring dda0feb
Merge branch 'main' into ExamplesNewtonKrylov
DanielDoehring 9497f15
Update test/Project.toml
DanielDoehring da5500b
Update Project.toml
DanielDoehring 73ee564
Merge branch 'main' into ExamplesNewtonKrylov
DanielDoehring 9820f59
Merge branch 'main' into ExamplesNewtonKrylov
DanielDoehring 8bbf603
Merge branch 'main' into ExamplesNewtonKrylov
DanielDoehring fcdfd39
ci
DanielDoehring 2365682
compat
DanielDoehring 6409d7c
v
DanielDoehring 98012be
v
DanielDoehring 178736e
v
DanielDoehring 1e8ae90
v
DanielDoehring 8a517fe
v
DanielDoehring bb9c185
v
DanielDoehring ae4d19d
Update test/Project.toml
DanielDoehring d7a4d77
Update test/Project.toml
DanielDoehring e04e649
Update test/Project.toml
JoshuaLampert add1a0c
Update test/Project.toml
JoshuaLampert e843d3b
Merge branch 'main' into ExamplesNewtonKrylov
DanielDoehring 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
171 changes: 171 additions & 0 deletions
171
examples/p4est_2d_dgsem/elixir_navierstokes_viscous_shock_newton_krylov.jl
This file contains hidden or 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 |
---|---|---|
@@ -0,0 +1,171 @@ | ||
using Trixi | ||
using OrdinaryDiffEqSDIRK | ||
using LinearSolve # For Jacobian-free Newton-Krylov (GMRES) solver | ||
using ADTypes # For automatic differentiation via finite differences | ||
|
||
# This is the classic 1D viscous shock wave problem with analytical solution | ||
# for a special value of the Prandtl number. | ||
# The original references are: | ||
# | ||
# - R. Becker (1922) | ||
# Stoßwelle und Detonation. | ||
# [DOI: 10.1007/BF01329605](https://doi.org/10.1007/BF01329605) | ||
# | ||
# English translations: | ||
# Impact waves and detonation. Part I. | ||
# https://ntrs.nasa.gov/api/citations/19930090862/downloads/19930090862.pdf | ||
# Impact waves and detonation. Part II. | ||
# https://ntrs.nasa.gov/api/citations/19930090863/downloads/19930090863.pdf | ||
# | ||
# - M. Morduchow, P. A. Libby (1949) | ||
# On a Complete Solution of the One-Dimensional Flow Equations | ||
# of a Viscous, Head-Conducting, Compressible Gas | ||
# [DOI: 10.2514/8.11882](https://doi.org/10.2514/8.11882) | ||
# | ||
# | ||
# The particular problem considered here is described in | ||
# - L. G. Margolin, J. M. Reisner, P. M. Jordan (2017) | ||
# Entropy in self-similar shock profiles | ||
# [DOI: 10.1016/j.ijnonlinmec.2017.07.003](https://doi.org/10.1016/j.ijnonlinmec.2017.07.003) | ||
|
||
### Fixed parameters ### | ||
|
||
# Special value for which nonlinear solver can be omitted | ||
# Corresponds essentially to fixing the Mach number | ||
alpha = 0.5 | ||
# We want kappa = cp * mu = mu_bar to ensure constant enthalpy | ||
prandtl_number() = 3 / 4 | ||
|
||
### Free choices: ### | ||
gamma() = 5 / 3 | ||
|
||
# In Margolin et al., the Navier-Stokes equations are given for an | ||
# isotropic stress tensor τ, i.e., ∇ ⋅ τ = μ Δu | ||
mu_isotropic() = 0.15 | ||
mu_bar() = mu_isotropic() / (gamma() - 1) # Re-scaled viscosity | ||
|
||
rho_0() = 1 | ||
v() = 1 # Shock speed | ||
|
||
domain_length = 4.0 | ||
|
||
### Derived quantities ### | ||
|
||
Ma() = 2 / sqrt(3 - gamma()) # Mach number for alpha = 0.5 | ||
c_0() = v() / Ma() # Speed of sound ahead of the shock | ||
|
||
# From constant enthalpy condition | ||
p_0() = c_0()^2 * rho_0() / gamma() | ||
|
||
l() = mu_bar() / (rho_0() * v()) * 2 * gamma() / (gamma() + 1) # Appropriate length scale | ||
|
||
""" | ||
initial_condition_viscous_shock(x, t, equations) | ||
Classic 1D viscous shock wave problem with analytical solution | ||
for a special value of the Prandtl number. | ||
The version implemented here is described in | ||
- L. G. Margolin, J. M. Reisner, P. M. Jordan (2017) | ||
Entropy in self-similar shock profiles | ||
[DOI: 10.1016/j.ijnonlinmec.2017.07.003](https://doi.org/10.1016/j.ijnonlinmec.2017.07.003) | ||
""" | ||
function initial_condition_viscous_shock(x, t, equations) | ||
y = x[1] - v() * t # Translated coordinate | ||
|
||
# Coordinate transformation. See eq. (33) in Margolin et al. (2017) | ||
chi = 2 * exp(y / (2 * l())) | ||
|
||
w = 1 + 1 / (2 * chi^2) * (1 - sqrt(1 + 2 * chi^2)) | ||
|
||
rho = rho_0() / w | ||
u = v() * (1 - w) | ||
p = p_0() * 1 / w * (1 + (gamma() - 1) / 2 * Ma()^2 * (1 - w^2)) | ||
|
||
return prim2cons(SVector(rho, u, 0, p), equations) | ||
end | ||
initial_condition = initial_condition_viscous_shock | ||
|
||
############################################################################### | ||
# semidiscretization of the ideal compressible Navier-Stokes equations | ||
|
||
equations = CompressibleEulerEquations2D(gamma()) | ||
|
||
# Trixi implements the stress tensor in deviatoric form, thus we need to | ||
# convert the "isotropic viscosity" to the "deviatoric viscosity" | ||
mu_deviatoric() = mu_bar() * 3 / 4 | ||
equations_parabolic = CompressibleNavierStokesDiffusion2D(equations, mu = mu_deviatoric(), | ||
Prandtl = prandtl_number(), | ||
gradient_variables = GradientVariablesEntropy()) | ||
|
||
solver = DGSEM(polydeg = 3, surface_flux = flux_hlle) | ||
|
||
coordinates_min = (-domain_length / 2, -domain_length / 2) | ||
coordinates_max = (domain_length / 2, domain_length / 2) | ||
|
||
trees_per_dimension = (12, 3) | ||
mesh = P4estMesh(trees_per_dimension, | ||
polydeg = 3, initial_refinement_level = 0, | ||
coordinates_min = coordinates_min, coordinates_max = coordinates_max, | ||
periodicity = (false, true)) | ||
|
||
### Inviscid boundary conditions ### | ||
|
||
# Prescribe pure influx based on initial conditions | ||
function boundary_condition_inflow(u_inner, normal_direction::AbstractVector, x, t, | ||
surface_flux_function, | ||
equations::CompressibleEulerEquations2D) | ||
u_cons = initial_condition_viscous_shock(x, t, equations) | ||
flux = Trixi.flux(u_cons, normal_direction, equations) | ||
|
||
return flux | ||
end | ||
|
||
boundary_conditions = Dict(:x_neg => boundary_condition_inflow, | ||
:x_pos => boundary_condition_do_nothing) | ||
|
||
### Viscous boundary conditions ### | ||
# For the viscous BCs, we use the known analytical solution | ||
velocity_bc = NoSlip() do x, t, equations_parabolic | ||
Trixi.velocity(initial_condition_viscous_shock(x, | ||
t, | ||
equations_parabolic), | ||
equations_parabolic) | ||
end | ||
|
||
heat_bc = Isothermal() do x, t, equations_parabolic | ||
Trixi.temperature(initial_condition_viscous_shock(x, | ||
t, | ||
equations_parabolic), | ||
equations_parabolic) | ||
end | ||
|
||
boundary_condition_parabolic = BoundaryConditionNavierStokesWall(velocity_bc, heat_bc) | ||
|
||
boundary_conditions_parabolic = Dict(:x_neg => boundary_condition_parabolic, | ||
:x_pos => boundary_condition_parabolic) | ||
|
||
semi = SemidiscretizationHyperbolicParabolic(mesh, (equations, equations_parabolic), | ||
initial_condition, solver; | ||
boundary_conditions = (boundary_conditions, | ||
boundary_conditions_parabolic)) | ||
|
||
############################################################################### | ||
# ODE solvers, callbacks etc. | ||
|
||
tspan = (0.0, 1.0) | ||
ode = semidiscretize(semi, tspan) | ||
|
||
summary_callback = SummaryCallback() | ||
alive_callback = AliveCallback(alive_interval = 1) | ||
analysis_callback = AnalysisCallback(semi, interval = 100) | ||
|
||
callbacks = CallbackSet(summary_callback, alive_callback, analysis_callback) | ||
|
||
############################################################################### | ||
# run the simulation | ||
|
||
sol = solve(ode, | ||
# Use (diagonally) implicit Runge-Kutta method with Jacobian-free (!) Newton-Krylov (GMRES) solver | ||
# See https://docs.sciml.ai/DiffEqDocs/stable/tutorials/advanced_ode_example/#Using-Jacobian-Free-Newton-Krylov | ||
KenCarp47(autodiff = AutoFiniteDiff(), linsolve = KrylovJL_GMRES()); | ||
ode_default_options()..., callback = callbacks) |
This file contains hidden or 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
59 changes: 59 additions & 0 deletions
59
examples/tree_1d_dgsem/elixir_diffusion_ldg_newton_krylov.jl
This file contains hidden or 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 |
---|---|---|
@@ -0,0 +1,59 @@ | ||
using Trixi | ||
using OrdinaryDiffEqSDIRK | ||
using LinearSolve # For Jacobian-free Newton-Krylov (GMRES) solver | ||
using ADTypes # For automatic differentiation via finite differences | ||
|
||
############################################################################### | ||
# semidiscretization of the linear (advection) diffusion equation | ||
|
||
advection_velocity = 0.0 # Note: This renders the equation mathematically purely parabolic | ||
equations = LinearScalarAdvectionEquation1D(advection_velocity) | ||
diffusivity() = 0.5 | ||
equations_parabolic = LaplaceDiffusion1D(diffusivity(), equations) | ||
|
||
# surface flux does not matter for pure diffusion problem | ||
solver = DGSEM(polydeg = 3, surface_flux = flux_central) | ||
|
||
coordinates_min = -convert(Float64, pi) | ||
coordinates_max = convert(Float64, pi) | ||
|
||
mesh = TreeMesh(coordinates_min, coordinates_max, | ||
initial_refinement_level = 4, | ||
n_cells_max = 30_000) | ||
|
||
function initial_condition_pure_diffusion_1d_convergence_test(x, t, | ||
equation) | ||
nu = diffusivity() | ||
c = 0 | ||
A = 1 | ||
omega = 1 | ||
scalar = c + A * sin(omega * sum(x)) * exp(-nu * omega^2 * t) | ||
return SVector(scalar) | ||
end | ||
initial_condition = initial_condition_pure_diffusion_1d_convergence_test | ||
|
||
solver_parabolic = ViscousFormulationLocalDG() | ||
semi = SemidiscretizationHyperbolicParabolic(mesh, (equations, equations_parabolic), | ||
initial_condition, | ||
solver; solver_parabolic) | ||
|
||
############################################################################### | ||
# ODE solvers, callbacks etc. | ||
|
||
tspan = (0.0, 2.0) | ||
ode = semidiscretize(semi, tspan) | ||
|
||
summary_callback = SummaryCallback() | ||
analysis_callback = AnalysisCallback(semi, interval = 10) | ||
alive_callback = AliveCallback(alive_interval = 1) | ||
|
||
callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback) | ||
|
||
############################################################################### | ||
# run the simulation | ||
|
||
sol = solve(ode, | ||
# Use (diagonally) implicit Runge-Kutta method with Jacobian-free (!) Newton-Krylov (GMRES) solver | ||
# See https://docs.sciml.ai/DiffEqDocs/stable/tutorials/advanced_ode_example/#Using-Jacobian-Free-Newton-Krylov | ||
KenCarp47(autodiff = AutoFiniteDiff(), linsolve = KrylovJL_GMRES()); | ||
ode_default_options()..., callback = callbacks) |
This file contains hidden or 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
This file contains hidden or 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
This file contains hidden or 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
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.
Uh oh!
There was an error while loading. Please reload this page.