-
Notifications
You must be signed in to change notification settings - Fork 121
Add example case for convergence test in 1D #1030
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
wilfonba
wants to merge
6
commits into
MFlowCode:master
Choose a base branch
from
wilfonba:convergenceExample
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.
Open
Changes from all commits
Commits
Show all changes
6 commits
Select commit
Hold shift + click to select a range
aad533c
add 1D convergence test example
wilfonba 87e0dcb
format
wilfonba f1c26f9
add missing files
wilfonba da82463
Merge branch 'master' into convergenceExample
wilfonba 3dd46e2
Github suggestions and bug fix
wilfonba 103ef16
Merge branch 'convergenceExample' of github.com:wilfonba/MFC-Wilfong …
wilfonba 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 |
|---|---|---|
| @@ -0,0 +1,11 @@ | ||
| This example case contains an automated convergence test using a 1D, two-component advection case. | ||
| The case can be run by executing the bash script `./submitJobs.sh` in a terminal after enabling | ||
| execution permissions with `chmod +x ./submitJobs.sh` and setting the `ROOT_DIR` and `MFC_DIR` | ||
| variables. By default the script runs the case for 6 different grid resolutions with 1st, | ||
| 3rd, and 5th, order spatial reconstructions. These settings can be modified by editing the variables | ||
| at the top of the script. You can also run different model equations by setting the `ME` variable | ||
| and different Riemann solvers by setting the `RS` variable. | ||
|
|
||
| Once the simulations have been run, you can generate convergence plots with matplotlib by running | ||
| `python3 plot.py` in a terminal. This will generate plots of the L1, L2, and Linf error norms and | ||
| save the results to `errors.csv`. |
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,81 @@ | ||
| #!/usr/bin/env python3 | ||
| import math | ||
| import json | ||
| import argparse | ||
|
|
||
| # Parsing command line arguments | ||
| parser = argparse.ArgumentParser(description="Generate JSON case configuration for two-fluid convergence simulation.") | ||
| parser.add_argument("--mfc", type=json.loads, default="{}", metavar="DICT", help="MFC's toolchain's internal state.") | ||
| parser.add_argument("--order", type=int, default=5, help="WENO order (default: 5)") | ||
| parser.add_argument("--meqns", type=int, default=2, help="Model equations (default: 2 (five-equation model))") | ||
| parser.add_argument("--rs", type=int, default=2, help="Riemann solver (default: 2 (HLLC))") | ||
| parser.add_argument("-N", type=int, default=1024, help="Number of grid points (default: 1024)") | ||
|
|
||
| args = parser.parse_args() | ||
|
|
||
| # Numerical setup | ||
| Nx = args.N - 1 | ||
| dx = 1.0 / (1.0 * (Nx + 1)) | ||
|
|
||
| Tend, Nt = 5.0, 200000 | ||
| mydt = Tend / (1.0 * Nt) | ||
|
|
||
| # Configuring case dictionary | ||
| print( | ||
| json.dumps( | ||
| { | ||
| # Logistics | ||
| "run_time_info": "T", | ||
| # Computational Domain Parameters | ||
| "x_domain%beg": 0.0e00, | ||
| "x_domain%end": 1.0e00, | ||
| "m": Nx, | ||
| "n": 0, | ||
| "p": 0, | ||
| "dt": mydt, | ||
| "t_step_start": 0, | ||
| "t_step_stop": int(Nt), | ||
| "t_step_save": int(math.ceil(Nt / 10.0)), | ||
| # Simulation Algorithm Parameters | ||
| "num_patches": 1, | ||
| "model_eqns": args.meqns, | ||
| "alt_soundspeed": "F", | ||
| "num_fluids": 2, | ||
| "mpp_lim": "F", | ||
| "mixture_err": "F", | ||
| "time_stepper": 3, | ||
| "weno_order": args.order, | ||
| "weno_eps": dx**2, | ||
| "weno_Re_flux": "F", | ||
| "weno_avg": "F", | ||
| "mapped_weno": "F", | ||
| "null_weights": "F", | ||
| "mp_weno": "F", | ||
| "riemann_solver": args.rs, | ||
| "wave_speeds": 1, | ||
| "avg_state": 2, | ||
| "bc_x%beg": -1, | ||
| "bc_x%end": -1, | ||
| # Formatted Database Files Structure Parameters | ||
| "format": 1, | ||
| "precision": 2, | ||
| "prim_vars_wrt": "T", | ||
| "parallel_io": "F", | ||
| # Patch 1 L | ||
| "patch_icpp(1)%geometry": 1, | ||
| "patch_icpp(1)%x_centroid": 0.5, | ||
| "patch_icpp(1)%length_x": 1.0, | ||
| "patch_icpp(1)%vel(1)": 1.0, | ||
| "patch_icpp(1)%pres": 1.0, | ||
| "patch_icpp(1)%alpha_rho(1)": f"0.5 - 0.5*sin(2*pi*x)", | ||
| "patch_icpp(1)%alpha(1)": f"0.5 - 0.5*sin(2*pi*x)", | ||
| "patch_icpp(1)%alpha_rho(2)": f"0.5 + 0.5*sin(2*pi*x)", | ||
| "patch_icpp(1)%alpha(2)": f"0.5 + 0.5*sin(2*pi*x)", | ||
| # Fluids Physical Parameters | ||
| "fluid_pp(1)%gamma": 1.0e00 / (1.4 - 1.0e00), | ||
| "fluid_pp(1)%pi_inf": 0.0, | ||
| "fluid_pp(2)%gamma": 1.0e00 / (1.4 - 1.0e00), | ||
| "fluid_pp(2)%pi_inf": 0.0, | ||
| } | ||
| ) | ||
| ) |
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,69 @@ | ||
| import matplotlib.pyplot as plt | ||
| import numpy as np | ||
| import pandas as pd | ||
|
|
||
| N = [32, 64, 128, 256, 512, 1024] | ||
| Ord = [1, 3, 5] | ||
|
|
||
| errors = np.nan * np.zeros((len(N), len(Ord), 3)) | ||
|
|
||
| TEND = 200000 | ||
|
|
||
| for i in range(len(N)): | ||
| for j in range(len(Ord)): | ||
|
|
||
| sim_a1 = pd.read_csv(f"N{N[i]}_O{Ord[j]}/D/cons.5.00.{TEND}.dat", sep=r"\s+", header=None, names=["x", "y"]) | ||
| sim_a2 = pd.read_csv(f"N{N[i]}_O{Ord[j]}/D/cons.6.00.{TEND}.dat", sep=r"\s+", header=None, names=["x", "y"]) | ||
|
|
||
| exact_a1 = pd.read_csv(f"N{N[i]}_O{Ord[j]}/D/cons.5.00.000000.dat", sep=r"\s+", header=None, names=["x", "y"]) | ||
| exact_a2 = pd.read_csv(f"N{N[i]}_O{Ord[j]}/D/cons.6.00.000000.dat", sep=r"\s+", header=None, names=["x", "y"]) | ||
|
|
||
| ## 2 norm | ||
| errors[i, j, 0] = np.linalg.norm(sim_a1.y - exact_a1.y) / np.sqrt(N[i]) | ||
| errors[i, j, 0] += np.linalg.norm(sim_a2.y - exact_a2.y) / np.sqrt(N[i]) | ||
|
|
||
| ## 1 norm | ||
| errors[i, j, 1] = 1 / N[i] * np.sum(np.abs(sim_a1.y - exact_a1.y)) | ||
| errors[i, j, 1] += 1 / N[i] * np.sum(np.abs(sim_a2.y - exact_a2.y)) | ||
|
|
||
| ## Inf norm | ||
| errors[i, j, 2] = np.max([np.nanmax(np.abs(sim_a1.y - exact_a1.y)), np.nanmax(np.abs(sim_a2.y - exact_a2.y))]) | ||
|
|
||
| fig, ax = plt.subplots(1, 3, figsize=(12, 8), sharex=True) | ||
|
|
||
| colors = ["blue", "green", "red", "purple"] | ||
|
|
||
| ref = np.nan * np.zeros((len(N), len(Ord))) | ||
|
|
||
| for i in range(3): | ||
| ax[i].plot(N, 30 / np.array(N) ** 1, label="Slope = -1", color=colors[0]) | ||
| ax[i].plot(N, 3000 / np.array(N) ** 3, label="Slope = -3", color=colors[1]) | ||
| ax[i].plot(N, 5000 / np.array(N) ** 5, label="Slope = -5", color=colors[2]) | ||
|
|
||
| for j in range(len(Ord)): | ||
| ax[0].plot(N, errors[:, j, 0], "o-", color=colors[j]) | ||
| ax[0].set_xscale("log", base=2) | ||
| ax[0].set_yscale("log") | ||
| ax[0].set_title("||error||_2") | ||
| ax[0].legend() | ||
|
|
||
| ax[1].plot(N, errors[:, j, 1], "o-", color=colors[j]) | ||
| ax[1].set_xscale("log", base=2) | ||
| ax[1].set_yscale("log") | ||
| ax[1].set_title("||error||_1") | ||
|
|
||
| ax[2].plot(N, errors[:, j, 2], "o-", color=colors[j]) | ||
| ax[2].set_xscale("log", base=2) | ||
| ax[2].set_yscale("log") | ||
| ax[2].set_title("||error||_inf") | ||
|
|
||
| plt.tight_layout() | ||
| plt.show() | ||
|
|
||
| errors = np.column_stack((N, errors[:, :, 0])) | ||
| errors = np.column_stack((errors, 7 / np.array(N) ** 1)) | ||
| errors = np.column_stack((errors, 700 / np.array(N) ** 3)) | ||
| errors = np.column_stack((errors, 2000 / np.array(N) ** 5)) | ||
|
|
||
| df = pd.DataFrame(errors, columns=["N", "1", "3", "5", "R1", "R3", "R5"], index=N) | ||
| df.to_csv("errors.csv", index=False) | ||
wilfonba marked this conversation as resolved.
Show resolved
Hide resolved
|
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,28 @@ | ||
| #!/usr/bin/env sh | ||
| Nx=(32 64 128 256 512 1024) | ||
| Order=(1 3 5) | ||
|
|
||
| ME=2 # Model equations = 2 for five-equation model | ||
| RS=2 # Riemann solver = 2 for HLLC | ||
|
|
||
| # Get the directory of the script itself | ||
| ROOT_DIR="$( cd "$( dirname "${BASH_SOURCE[0]}" )" &> /dev/null && pwd )" | ||
| # Assume the script is in examples/1D_convergence, so MFC_DIR is two levels up | ||
| MFC_DIR="$(dirname "$(dirname "$ROOT_DIR")")" | ||
|
|
||
| for i in "${Nx[@]}"; do | ||
| for j in "${Order[@]}"; do | ||
| rm -rf N${i}_O${j} | ||
| mkdir N${i}_O${j} | ||
| cp case.py N${i}_O${j}/ | ||
| done | ||
| done | ||
|
|
||
| cd $MFC_DIR | ||
|
|
||
| for i in "${Nx[@]}"; do | ||
| for j in "${Order[@]}"; do | ||
| ./mfc.sh run $ROOT_DIR/N${i}_O${j}/case.py --case-optimization --no-debug -- --order $j -N $i --meqns $ME --rs $RS | ||
| done | ||
| done | ||
|
|
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.