-
-
Notifications
You must be signed in to change notification settings - Fork 128
Add julia example in resonant circuit #658
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
marinlauber
wants to merge
5
commits into
precice:develop
Choose a base branch
from
marinlauber:add-julia-resonant-circuit
base: develop
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 all commits
Commits
Show all changes
5 commits
Select commit
Hold shift + click to select a range
a5f9532
add julia resonant circuit example
marinlauber 114d811
update README.md
marinlauber 64889d9
tidy scripts and remove plotting from there, fix run.sh
marinlauber 6d4bae9
add package install in run and switch OrdinaryDiffEq for Differential…
marinlauber cc4dedf
rm old Project files
marinlauber 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 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 |
---|---|---|
@@ -1,7 +1,7 @@ | ||
--- | ||
title: Resonant Circuit | ||
permalink: tutorials-resonant-circuit.html | ||
keywords: MATLAB, Python | ||
keywords: MATLAB, Python, Julia | ||
summary: We simulate a two-element LC circuit (one inductor and one capacitor). | ||
--- | ||
|
||
|
@@ -33,6 +33,7 @@ preCICE configuration (image generated using the [precice-config-visualizer](htt | |
* *MATLAB* A solver using the [MATLAB bindings](https://precice.org/installation-bindings-matlab.html). | ||
Before running this tutorial, follow the [instructions](https://precice.org/installation-bindings-matlab.html) to correctly install the MATLAB bindings. | ||
* *Python* A solver using the preCICE [Python bindings](https://precice.org/installation-bindings-python.html). | ||
* *Julia* A solver using the preCICE [Julia bindings](https://precice.org/installation-bindings-julia.html). | ||
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. Self-reminder: We will need also an edit in https://precice.org/tutorials.html |
||
|
||
## Running the simulation | ||
|
||
|
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,82 @@ | ||
using PreCICE | ||
using DifferentialEquations | ||
|
||
# Initialize and configure preCICE | ||
participant = PreCICE.createParticipant("Capacitor", "../precice-config.xml", 0, 1) | ||
|
||
# Geometry IDs. As it is a 0-D simulation, only one vertex is necessary. | ||
mesh_name = "Capacitor-Mesh" | ||
|
||
dimensions = PreCICE.getMeshDimensions(mesh_name) | ||
|
||
vertex_ids = PreCICE.setMeshVertices(mesh_name, zeros((1,dimensions))) | ||
|
||
let | ||
# Data IDs | ||
read_data_name = "Current" | ||
write_data_name = "Voltage" | ||
|
||
# Simulation parameters and initial condition | ||
C = 2 # Capacitance of the capacitor | ||
L = 1 # Inductance of the coil | ||
t0 = 0 # Initial simulation time | ||
t_max = 10 # End simulation time | ||
Io = 1 # Initial current | ||
phi = 0 # Phase of the signal | ||
|
||
w0 = 1 / sqrt(L * C) # Resonant frequency | ||
U0 = -w0 * L * Io * sin(phi) # Initial condition for U | ||
|
||
function f_U(u, p, t) | ||
(dt, C, mesh_name, read_data_name, vertex_ids) = p | ||
I = PreCICE.readData(mesh_name, read_data_name, vertex_ids, dt) | ||
return -I[1] / C # Time derivative of U | ||
end | ||
MakisH marked this conversation as resolved.
Show resolved
Hide resolved
|
||
|
||
# Initialize simulation | ||
if PreCICE.requiresInitialData() | ||
PreCICE.writeData(mesh_name, write_data_name, vertex_ids, I0) | ||
end | ||
PreCICE.initialize() | ||
|
||
solver_dt = PreCICE.getMaxTimeStepSize() | ||
|
||
# Start simulation | ||
t = t0 | ||
U0_checkpoint = U0 | ||
t_checkpoint = t | ||
while PreCICE.isCouplingOngoing() | ||
|
||
# Record checkpoint if necessary | ||
if PreCICE.requiresWritingCheckpoint() | ||
U0_checkpoint = U0 | ||
t_checkpoint = t | ||
end | ||
|
||
# Make Simulation Step | ||
precice_dt = PreCICE.getMaxTimeStepSize() | ||
dt = min(precice_dt, solver_dt) | ||
t_span = (t, t + dt) | ||
params = (dt, C, mesh_name, read_data_name, vertex_ids) | ||
prob = ODEProblem(f_U, U0, t_span, params) | ||
# https://docs.sciml.ai/DiffEqDocs/stable/basics/common_solver_opts | ||
sol = solve(prob, Tsit5(), reltol=1e-12, abstol=1e-12) | ||
|
||
# Exchange data | ||
evals = max(length(sol.t), 3) # at least do 3 substeps to allow cubic interpolation | ||
for i in range(1,evals) | ||
U0 = sol(t_span[1] + i * dt / evals) | ||
PreCICE.writeData(mesh_name, write_data_name, vertex_ids, fill(U0, (1,1))) | ||
PreCICE.advance(dt / evals) | ||
end | ||
t = t + dt | ||
|
||
# Recover checkpoint if not converged | ||
if PreCICE.requiresReadingCheckpoint() | ||
U0 = U0_checkpoint | ||
t = t_checkpoint | ||
end | ||
end | ||
# Stop coupling | ||
PreCICE.finalize() | ||
end |
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,10 @@ | ||
#!/usr/bin/env bash | ||
set -e -u | ||
|
||
. ../../tools/log.sh | ||
exec > >(tee --append "$LOGFILE") 2>&1 | ||
|
||
julia --project=Project.toml -e 'using Pkg; Pkg.add("DifferentialEquations"); Pkg.add("PreCICE"); Pkg.instantiate();' | ||
julia --project=Project.toml capacitor.jl | ||
|
||
close_log |
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,85 @@ | ||
using PreCICE | ||
using DifferentialEquations | ||
|
||
# Initialize and configure preCICE | ||
participant = PreCICE.createParticipant("Coil", "../precice-config.xml", 0, 1) | ||
|
||
# Geometry IDs. As it is a 0-D simulation, only one vertex is necessary. | ||
mesh_name = "Coil-Mesh" | ||
|
||
dimensions = PreCICE.getMeshDimensions(mesh_name) | ||
|
||
vertex_ids = PreCICE.setMeshVertices(mesh_name, zeros((1,dimensions))) | ||
|
||
let | ||
# Data IDs | ||
read_data_name = "Voltage" | ||
write_data_name = "Current" | ||
|
||
# Simulation parameters and initial condition | ||
C = 2 # Capacitance of the capacitor | ||
L = 1 # Inductance of the coil | ||
t0 = 0 # Initial simulation time | ||
Io = 1 # Initial current | ||
phi = 0 # Phase of the signal | ||
|
||
w0 = 1 / sqrt(L * C) # Resonant frequency | ||
I0 = Io * cos(phi) # Initial condition for I | ||
|
||
# to estimate cost | ||
f_evals = 0 | ||
|
||
function f_I(u, p, t) | ||
f_evals += 1 | ||
(dt, L, mesh_name, read_data_name, vertex_ids) = p | ||
U = PreCICE.readData(mesh_name, read_data_name, vertex_ids, dt) | ||
return U[1] / L # Time derivative of I; ODE determining capacitor | ||
end | ||
|
||
# Initialize simulation | ||
if PreCICE.requiresInitialData() | ||
PreCICE.writeData(mesh_name, write_data_name, vertex_ids, I0) | ||
end | ||
PreCICE.initialize() | ||
|
||
solver_dt = PreCICE.getMaxTimeStepSize() | ||
|
||
# Start simulation | ||
t = t0 | ||
I0_checkpoint = I0 | ||
t_checkpoint = t | ||
while PreCICE.isCouplingOngoing() | ||
|
||
# Record checkpoint if necessary | ||
if PreCICE.requiresWritingCheckpoint() | ||
I0_checkpoint = I0 | ||
t_checkpoint = t | ||
end | ||
|
||
# Make Simulation Step | ||
precice_dt = PreCICE.getMaxTimeStepSize() | ||
dt = min(precice_dt, solver_dt) | ||
t_span = (t, t + dt) | ||
params = (dt, L, mesh_name, read_data_name, vertex_ids) | ||
prob = ODEProblem(f_I, I0, t_span, params) | ||
# https://docs.sciml.ai/DiffEqDocs/stable/basics/common_solver_opts | ||
sol = solve(prob, Tsit5(), reltol=1e-12, abstol=1e-12) | ||
|
||
# Exchange data | ||
evals = max(length(sol.t), 3) # at least do 3 substeps to allow cubic interpolation | ||
for i in range(1,evals) | ||
I0 = sol(t_span[1] + i * dt / evals) | ||
PreCICE.writeData(mesh_name, write_data_name, vertex_ids, fill(I0, (1,1))) | ||
PreCICE.advance(dt / evals) | ||
end | ||
t = t + dt | ||
|
||
# Recover checkpoint if not converged | ||
if PreCICE.requiresReadingCheckpoint() | ||
I0 = I0_checkpoint | ||
t = t_checkpoint | ||
end | ||
end | ||
# Stop coupling | ||
PreCICE.finalize() | ||
end |
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,10 @@ | ||
#!/usr/bin/env bash | ||
set -e -u | ||
|
||
. ../../tools/log.sh | ||
exec > >(tee --append "$LOGFILE") 2>&1 | ||
|
||
julia --project=Project.toml -e 'using Pkg; Pkg.add("DifferentialEquations"); Pkg.add("PreCICE"); Pkg.instantiate();' | ||
julia --project=Project.toml coil.jl | ||
|
||
close_log |
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.
The rest of the page seems to be a bit focused on MATLAB. This is a good chance to change this (I can do).