From 57bf01e4db3f006da866d0a0f6df7561b3175933 Mon Sep 17 00:00:00 2001 From: Dennis Yatunin Date: Tue, 10 Jan 2023 11:23:59 -0800 Subject: [PATCH] Add integration test with limiters --- .buildkite/pipeline.yml | 2 - docs/Manifest.toml | 115 ++++++++++---- docs/Project.toml | 2 + docs/src/dev/report_gen.jl | 28 ++-- docs/src/dev/report_gen.md | 43 +++++- docs/src/plotting_utils.jl | 131 +++++++++++++--- src/solvers/imex_ssp.jl | 4 +- test/problems.jl | 302 +++++++++++++++++++++++++++++++++++++ 8 files changed, 555 insertions(+), 72 deletions(-) diff --git a/.buildkite/pipeline.yml b/.buildkite/pipeline.yml index 2525bdd9..5aa39462 100644 --- a/.buildkite/pipeline.yml +++ b/.buildkite/pipeline.yml @@ -4,13 +4,11 @@ env: OPENMPI_VERSION: "4.1.1" OPENBLAS_NUM_THREADS: 1 CLIMATEMACHINE_SETTINGS_FIX_RNG_SEED: "true" - JULIA_DEPOT_PATH: "${BUILDKITE_BUILD_PATH}/${BUILDKITE_PIPELINE_SLUG}/depot/cpu" steps: - label: "init cpu env" key: "init_cpu_env" command: - - "echo $$JULIA_DEPOT_PATH" - echo "--- Configure MPI" - julia -e 'using Pkg; Pkg.add("MPIPreferences"); using MPIPreferences; use_system_binary()' diff --git a/docs/Manifest.toml b/docs/Manifest.toml index f60e2f22..f6eb6c8a 100644 --- a/docs/Manifest.toml +++ b/docs/Manifest.toml @@ -2,7 +2,7 @@ julia_version = "1.8.5" manifest_format = "2.0" -project_hash = "b41a2b6962844aa850e2210be1742cc1868a096b" +project_hash = "0445f0b63f403d783217e170583dfd95f7a70fe5" [[deps.AMD]] deps = ["Libdl", "LinearAlgebra", "SparseArrays", "Test"] @@ -173,9 +173,9 @@ version = "0.5.1" [[deps.ChainRulesCore]] deps = ["Compat", "LinearAlgebra", "SparseArrays"] -git-tree-sha1 = "e7ff6cadf743c098e08fca25c91103ee4303c9bb" +git-tree-sha1 = "c6d890a52d2c4d55d326439580c3b8d0875a77d9" uuid = "d360d2e6-b24c-11e9-a2a3-2a2ae2dbcce4" -version = "1.15.6" +version = "1.15.7" [[deps.ChangesOfVariables]] deps = ["ChainRulesCore", "LinearAlgebra", "Test"] @@ -195,6 +195,12 @@ git-tree-sha1 = "ba403ef47480ab31edb6b02d499538574b321c74" uuid = "d414da3d-4745-48bb-8d80-42e94e092884" version = "0.10.23" +[[deps.ClimaCorePlots]] +deps = ["ClimaCore", "RecipesBase", "StaticArrays", "TriplotBase"] +git-tree-sha1 = "2fb63d2573c4ba4aff1ed7829a6c91032be5d488" +uuid = "cf7c7e5a-b407-4c48-9047-11a94a308626" +version = "0.2.4" + [[deps.ClimaTimeSteppers]] deps = ["CUDA", "ClimaComms", "DataStructures", "DiffEqBase", "DiffEqCallbacks", "KernelAbstractions", "Krylov", "LinearAlgebra", "LinearOperators", "SciMLBase", "StaticArrays"] path = ".." @@ -276,11 +282,16 @@ git-tree-sha1 = "fcbb72b032692610bfbdb15018ac16a36cf2e406" uuid = "adafc99b-e345-5852-983c-f28acb93d879" version = "0.3.1" +[[deps.Crayons]] +git-tree-sha1 = "249fe38abf76d48563e2f4556bebd215aa317e15" +uuid = "a8cc5b0e-0ffa-5ad4-8c14-923d3ee1735f" +version = "4.1.1" + [[deps.CubedSphere]] -deps = ["Elliptic", "Printf", "Rotations", "TaylorSeries", "Test"] -git-tree-sha1 = "d5b15ff3bc2073a864f38a09b446a3b735a96125" +deps = ["Elliptic", "FFTW", "Printf", "ProgressBars", "SpecialFunctions", "TaylorSeries", "Test"] +git-tree-sha1 = "db9c12cb765cc048e158987388287c52baddf57d" uuid = "7445602f-e544-4518-8976-18f8e8ae6cdb" -version = "0.2.1" +version = "0.2.2" [[deps.DataAPI]] git-tree-sha1 = "e8119c1a33d267e16108be441a287a6981ba1630" @@ -420,6 +431,18 @@ git-tree-sha1 = "74faea50c1d007c85837327f6775bea60b5492dd" uuid = "b22a6f82-2f65-5046-a5b2-351ab43fb4e5" version = "4.4.2+2" +[[deps.FFTW]] +deps = ["AbstractFFTs", "FFTW_jll", "LinearAlgebra", "MKL_jll", "Preferences", "Reexport"] +git-tree-sha1 = "90630efff0894f8142308e334473eba54c433549" +uuid = "7a1cc6ca-52ef-59f5-83cd-3a7055c09341" +version = "1.5.0" + +[[deps.FFTW_jll]] +deps = ["Artifacts", "JLLWrappers", "Libdl", "Pkg"] +git-tree-sha1 = "c6033cc3892d0ef5bb9cd29b7f2f0331ea5184ea" +uuid = "f5851436-0d7a-5f13-b9de-f02708fd171a" +version = "3.3.10+0" + [[deps.FastBroadcast]] deps = ["ArrayInterface", "ArrayInterfaceCore", "LinearAlgebra", "Polyester", "Static", "StrideArraysCore"] git-tree-sha1 = "4bef892787c972913d4d84e7255400759bb650e5" @@ -517,15 +540,15 @@ version = "3.3.8+0" [[deps.GPUArrays]] deps = ["Adapt", "GPUArraysCore", "LLVM", "LinearAlgebra", "Printf", "Random", "Reexport", "Serialization", "Statistics"] -git-tree-sha1 = "45d7deaf05cbb44116ba785d147c518ab46352d7" +git-tree-sha1 = "5ae0f2ab6f6d99d558cb588763f76311106c71b1" uuid = "0c68f7d7-f131-5f86-a1c3-88cf8149b2d7" -version = "8.5.0" +version = "8.6.0" [[deps.GPUArraysCore]] deps = ["Adapt"] -git-tree-sha1 = "6872f5ec8fd1a38880f027a26739d42dcda6691f" +git-tree-sha1 = "57f7cde02d7a53c9d1d28443b9f11ac5fbe7ebc9" uuid = "46192b85-c4d5-4398-a991-12ede77f4527" -version = "0.1.2" +version = "0.1.3" [[deps.GPUCompiler]] deps = ["ExprTools", "InteractiveUtils", "LLVM", "Libdl", "Logging", "TimerOutputs", "UUIDs"] @@ -605,9 +628,9 @@ version = "1.12.2+2" [[deps.HTTP]] deps = ["Base64", "CodecZlib", "Dates", "IniFile", "Logging", "LoggingExtras", "MbedTLS", "NetworkOptions", "OpenSSL", "Random", "SimpleBufferStream", "Sockets", "URIs", "UUIDs"] -git-tree-sha1 = "fd9861adba6b9ae4b42582032d0936d456c8602d" +git-tree-sha1 = "eb5aa5e3b500e191763d35198f859e4b40fff4a6" uuid = "cd3eb016-35fb-5094-929b-558a96fad6f3" -version = "1.6.3" +version = "1.7.3" [[deps.HarfBuzz_jll]] deps = ["Artifacts", "Cairo_jll", "Fontconfig_jll", "FreeType2_jll", "Glib_jll", "Graphite2_jll", "JLLWrappers", "Libdl", "Libffi_jll", "Pkg"] @@ -648,6 +671,12 @@ git-tree-sha1 = "f550e6e32074c939295eb5ea6de31849ac2c9625" uuid = "83e8ac13-25f8-5344-8a64-a9f2b223428f" version = "0.5.1" +[[deps.IntelOpenMP_jll]] +deps = ["Artifacts", "JLLWrappers", "Libdl", "Pkg"] +git-tree-sha1 = "d979e54b71da82f3a65b62553da4fc3d18c9004c" +uuid = "1d5cc7b8-4909-519e-a0f8-d0f5ad9712d0" +version = "2018.0.3+2" + [[deps.InteractiveUtils]] deps = ["Markdown"] uuid = "b77e0a4c-d291-57a0-90e8-8db25a27a240" @@ -754,9 +783,9 @@ version = "3.0.0+1" [[deps.LLVM]] deps = ["CEnum", "LLVMExtra_jll", "Libdl", "Printf", "Unicode"] -git-tree-sha1 = "088dd02b2797f0233d92583562ab669de8517fd1" +git-tree-sha1 = "b8ae281340f0d3e973aae7b96fb7502b0119b376" uuid = "929cbde3-209d-540e-8aea-75f648917ca0" -version = "4.14.1" +version = "4.15.0" [[deps.LLVMExtra_jll]] deps = ["Artifacts", "JLLWrappers", "LazyArtifacts", "Libdl", "Pkg", "TOML"] @@ -900,9 +929,15 @@ version = "1.0.0" [[deps.LoopVectorization]] deps = ["ArrayInterface", "ArrayInterfaceCore", "ArrayInterfaceOffsetArrays", "ArrayInterfaceStaticArrays", "CPUSummary", "ChainRulesCore", "CloseOpenIntervals", "DocStringExtensions", "ForwardDiff", "HostCPUFeatures", "IfElse", "LayoutPointers", "LinearAlgebra", "OffsetArrays", "PolyesterWeave", "SIMDDualNumbers", "SIMDTypes", "SLEEFPirates", "SnoopPrecompile", "SpecialFunctions", "Static", "ThreadingUtilities", "UnPack", "VectorizationBase"] -git-tree-sha1 = "155132d68bc33c826dbdeb452c5d0a79e2d0e586" +git-tree-sha1 = "641ba2dbd7667d1ede0e9135aa3585018581d99b" uuid = "bdcacae8-1622-11e9-2a5c-532679323890" -version = "0.12.146" +version = "0.12.147" + +[[deps.MKL_jll]] +deps = ["Artifacts", "IntelOpenMP_jll", "JLLWrappers", "LazyArtifacts", "Libdl", "Pkg"] +git-tree-sha1 = "2ce8695e1e699b68702c03402672a69f54b8aca9" +uuid = "856f044c-d86e-5d09-b602-aeab76dc8ba7" +version = "2022.2.0+0" [[deps.MacroTools]] deps = ["Markdown", "Random"] @@ -1005,9 +1040,9 @@ version = "0.8.1+0" [[deps.OpenSSL]] deps = ["BitFlags", "Dates", "MozillaCACerts_jll", "OpenSSL_jll", "Sockets"] -git-tree-sha1 = "df6830e37943c7aaa10023471ca47fb3065cc3c4" +git-tree-sha1 = "6503b77492fd7fcb9379bf73cd31035670e3c509" uuid = "4d8831e6-92b7-49fb-bdf8-b643e874388c" -version = "1.3.2" +version = "1.3.3" [[deps.OpenSSL_jll]] deps = ["Artifacts", "JLLWrappers", "Libdl", "Pkg"] @@ -1057,9 +1092,9 @@ version = "0.12.3" [[deps.Parsers]] deps = ["Dates", "SnoopPrecompile"] -git-tree-sha1 = "6466e524967496866901a78fca3f2e9ea445a559" +git-tree-sha1 = "8175fc2b118a3755113c8e68084dc1a9e63c61ee" uuid = "69de0a69-1ddd-5017-9359-2bf0b02dc9f0" -version = "2.5.2" +version = "2.5.3" [[deps.Pipe]] git-tree-sha1 = "6842804e7867b115ca9de748a0cf6b364523c16d" @@ -1079,9 +1114,9 @@ version = "1.8.0" [[deps.PkgVersion]] deps = ["Pkg"] -git-tree-sha1 = "f6cf8e7944e50901594838951729a1861e668cb8" +git-tree-sha1 = "192aa02cc47f6566c7a1a8907dee76f6eb0a7c90" uuid = "eebad327-c553-4316-9ea0-9fa01ccd7688" -version = "0.3.2" +version = "0.2.1" [[deps.PlotThemes]] deps = ["PlotUtils", "Statistics"] @@ -1125,10 +1160,22 @@ git-tree-sha1 = "47e5f437cc0e7ef2ce8406ce1e7e24d44915f88d" uuid = "21216c6a-2e73-6563-6e65-726566657250" version = "1.3.0" +[[deps.PrettyTables]] +deps = ["Crayons", "Formatting", "LaTeXStrings", "Markdown", "Reexport", "StringManipulation", "Tables"] +git-tree-sha1 = "96f6db03ab535bdb901300f88335257b0018689d" +uuid = "08abe8d2-0d0c-5749-adfa-8a2ac140af0d" +version = "2.2.2" + [[deps.Printf]] deps = ["Unicode"] uuid = "de0858da-6303-5e67-8744-51eddeeeb8d7" +[[deps.ProgressBars]] +deps = ["Printf"] +git-tree-sha1 = "806ebc92e1b4b4f72192369a28dfcaf688566b2b" +uuid = "49802e3a-d2f1-5c88-81d8-b72133a6f568" +version = "1.4.1" + [[deps.Qt5Base_jll]] deps = ["Artifacts", "CompilerSupportLibraries_jll", "Fontconfig_jll", "Glib_jll", "JLLWrappers", "Libdl", "Libglvnd_jll", "OpenSSL_jll", "Pkg", "Xorg_libXext_jll", "Xorg_libxcb_jll", "Xorg_xcb_util_image_jll", "Xorg_xcb_util_keysyms_jll", "Xorg_xcb_util_renderutil_jll", "Xorg_xcb_util_wm_jll", "Zlib_jll", "xkbcommon_jll"] git-tree-sha1 = "0c03844e2231e12fda4d0086fd7cbe4098ee8dc5" @@ -1137,9 +1184,9 @@ version = "5.15.3+2" [[deps.QuadGK]] deps = ["DataStructures", "LinearAlgebra"] -git-tree-sha1 = "97aa253e65b784fd13e83774cadc95b38011d734" +git-tree-sha1 = "de191bc385072cc6c7ed3ffdc1caeed3f22c74d4" uuid = "1fd47b50-473d-5c70-9696-f719f8f3bcdc" -version = "2.6.0" +version = "2.7.0" [[deps.Quaternions]] deps = ["LinearAlgebra", "Random", "RealDot"] @@ -1193,9 +1240,9 @@ version = "2.35.0" [[deps.RecursiveFactorization]] deps = ["LinearAlgebra", "LoopVectorization", "Polyester", "SnoopPrecompile", "StrideArraysCore", "TriangularSolve"] -git-tree-sha1 = "9f9d83b485b5f3bfd0f8cb0b8733d573bd4c388f" +git-tree-sha1 = "664aba1c5259821356f2ef771eabc502d67a8f0d" uuid = "f2c3362d-daeb-58d1-803e-2bc74f2840b4" -version = "0.2.14" +version = "0.2.16" [[deps.Reexport]] git-tree-sha1 = "45e428421666073eab6f2da5c9d310d99bb12f9b" @@ -1308,10 +1355,10 @@ uuid = "777ac1f9-54b0-4bf8-805c-2214025038e7" version = "1.1.0" [[deps.SimpleNonlinearSolve]] -deps = ["ArrayInterfaceCore", "FiniteDiff", "ForwardDiff", "Reexport", "SciMLBase", "SnoopPrecompile", "StaticArraysCore"] -git-tree-sha1 = "fc4b9f81a033cf6879c91bb7f5b3ff59008c7dd2" +deps = ["ArrayInterfaceCore", "FiniteDiff", "ForwardDiff", "LinearAlgebra", "Reexport", "SciMLBase", "SnoopPrecompile", "StaticArraysCore"] +git-tree-sha1 = "61b8ffdb22453132e02a10c5638dfb42834c776b" uuid = "727e6d20-b764-4bd8-a329-72de5adea6c7" -version = "0.1.4" +version = "0.1.5" [[deps.SimpleTraits]] deps = ["InteractiveUtils", "MacroTools"] @@ -1407,6 +1454,11 @@ git-tree-sha1 = "50ccd5ddb00d19392577902f0079267a72c5ab04" uuid = "69024149-9ee7-55f6-a4c4-859efe599b68" version = "0.3.5" +[[deps.StringManipulation]] +git-tree-sha1 = "46da2434b41f41ac3594ee9816ce5541c6096123" +uuid = "892a3eda-7b42-436c-8928-eab12a02cf0e" +version = "0.3.0" + [[deps.SuiteSparse]] deps = ["Libdl", "LinearAlgebra", "Serialization", "SparseArrays"] uuid = "4607b0f0-06f3-5cda-b6b1-a6196a1729e9" @@ -1489,6 +1541,11 @@ git-tree-sha1 = "6bac775f2d42a611cdfcd1fb217ee719630c4175" uuid = "410a4b4d-49e4-4fbc-ab6d-cb71b17b3775" version = "0.1.6" +[[deps.TriplotBase]] +git-tree-sha1 = "4d4ed7f294cda19382ff7de4c137d24d16adc89b" +uuid = "981d1d27-644d-49a2-9326-4793e63143c3" +version = "0.1.0" + [[deps.URIs]] git-tree-sha1 = "ac00576f90d8a259f2c9d823e91d1de3fd44d348" uuid = "5c2747f8-b7ea-4ff2-ba2e-563bfd36b1d4" diff --git a/docs/Project.toml b/docs/Project.toml index b86352ab..63d69088 100644 --- a/docs/Project.toml +++ b/docs/Project.toml @@ -1,6 +1,7 @@ [deps] AbstractTrees = "1520ce14-60c1-5f80-bbc7-55ef81b5835c" ClimaCore = "d414da3d-4745-48bb-8d80-42e94e092884" +ClimaCorePlots = "cf7c7e5a-b407-4c48-9047-11a94a308626" ClimaTimeSteppers = "595c0a79-7f3d-439a-bc5a-b232dc3bde79" DiffEqBase = "2b5f629d-d688-5b77-993f-72d75c75574e" Distributions = "31c24e10-a181-5473-b8eb-7969acd0382f" @@ -9,5 +10,6 @@ DocumenterCitations = "daee34ce-89f3-4625-b898-19384cb65244" LaTeXStrings = "b964fa9f-0449-5b57-a5c2-d3ea65f4040f" OrdinaryDiffEq = "1dea7af3-3e70-54e6-95c3-0bf5283fa5ed" Plots = "91a5bcdd-55d7-5caf-9e0b-520d859cae80" +PrettyTables = "08abe8d2-0d0c-5749-adfa-8a2ac140af0d" Printf = "de0858da-6303-5e67-8744-51eddeeeb8d7" StaticArrays = "90137ffa-7385-5640-81b9-e52037218182" diff --git a/docs/src/dev/report_gen.jl b/docs/src/dev/report_gen.jl index ebdecc17..2226a99e 100644 --- a/docs/src/dev/report_gen.jl +++ b/docs/src/dev/report_gen.jl @@ -9,20 +9,14 @@ include(joinpath(pkgdir(ClimaTimeSteppers), "test", "problems.jl")) all_subtypes(::Type{T}) where {T} = isabstracttype(T) ? vcat(all_subtypes.(subtypes(T))...) : [T] -@testset "IMEX Algorithm Convergence" begin +let # IMEX Algorithm convergence title = "IMEX Algorithms" algorithm_names = map(T -> T(), all_subtypes(ClimaTimeSteppers.IMEXAlgorithmName)) - test_imex_algorithms(title, algorithm_names, ark_analytic_nonlin_test_cts(Float64), 200) - test_imex_algorithms(title, algorithm_names, ark_analytic_sys_test_cts(Float64), 400) - test_imex_algorithms(title, algorithm_names, ark_analytic_test_cts(Float64), 40000; super_convergence = (ARS121(),)) - test_imex_algorithms( - title, - algorithm_names, - onewaycouple_mri_test_cts(Float64), - 10000; - num_steps_scaling_factor = 5, - ) - test_imex_algorithms( + verify_convergence(title, algorithm_names, ark_analytic_nonlin_test_cts(Float64), 200) + verify_convergence(title, algorithm_names, ark_analytic_sys_test_cts(Float64), 400) + verify_convergence(title, algorithm_names, ark_analytic_test_cts(Float64), 40000; super_convergence = (ARS121(),)) + verify_convergence(title, algorithm_names, onewaycouple_mri_test_cts(Float64), 10000; num_steps_scaling_factor = 5) + verify_convergence( title, algorithm_names, climacore_1Dheat_test_cts(Float64), @@ -30,7 +24,7 @@ all_subtypes(::Type{T}) where {T} = isabstracttype(T) ? vcat(all_subtypes.(subty num_steps_scaling_factor = 4, numerical_reference_algorithm_name = ARS343(), ) - test_imex_algorithms( + verify_convergence( title, algorithm_names, climacore_2Dheat_test_cts(Float64), @@ -40,7 +34,7 @@ all_subtypes(::Type{T}) where {T} = isabstracttype(T) ? vcat(all_subtypes.(subty ) end -@testset "Unconstrained vs SSPConstrained (no limiters)" begin +let # Unconstrained vs SSPConstrained results without limiters algorithm_names = map(T -> T(), all_subtypes(ClimaTimeSteppers.IMEXSSPRKAlgorithmName)) for (test_case, num_steps) in ( (ark_analytic_nonlin_test_cts(Float64), 200), @@ -58,7 +52,11 @@ end reference_algorithm = IMEXAlgorithm(algorithm_name, newtons_method, Unconstrained()) solution = solve(deepcopy(prob), algorithm; dt).u[end] reference_solution = solve(deepcopy(prob), reference_algorithm; dt).u[end] - @test norm(solution .- reference_solution) / norm(reference_solution) < 30 * eps(Float64) + if norm(solution .- reference_solution) / norm(reference_solution) > 30 * eps(Float64) + alg_str = string(nameof(typeof(algorithm_name))) + @warn "Unconstrained and SSPConstrained versions of $alg_str \ + give different results for $(test_case.test_name)" + end end end end diff --git a/docs/src/dev/report_gen.md b/docs/src/dev/report_gen.md index 30835e83..353af379 100644 --- a/docs/src/dev/report_gen.md +++ b/docs/src/dev/report_gen.md @@ -1,19 +1,50 @@ -# Report generator +# Verifying Correctness -In this section, we create a report comparing the solution errors and convergence orders for several test problems and algorithms. +The `IMEXAlgorithm` supports problems that specify any combination of the following: an implicit tendency `T_imp!`, an explicit tendency `T_exp!`, a limited tendency `T_lim!`, a function `dss!` that applies a direct stiffness summation, and a function `lim!` that applies a monotonicity-preserving limiter. + +## Convergence without a Limiter + +In order to verify the correctness of our algorithms without a limiter, we compute their convergence orders for a variety of test cases. For each case, we estimate the convergence order of the algorithm over a range of stable timesteps, ensuring that the estimate `computed_order ± order_uncertainty` satisfies `|computed_order - predicted_order| ≤ order_uncertainty` and `order_uncertainty ≤ predicted_order / 10`. We also generate a plot that shows each algorithm's convergence as the timestep is reduced, along with plots that show the norms of each algorithm's solution and error over time (for some stable timestep). In addition, we verify that `SSPConstrained` algorithms produce the same results (up to floating-point roundoff error) when run in `Unconstrained` mode (at least, when run without a limiter). + +By [Godunov's theorem](https://en.wikipedia.org/wiki/Godunov%27s_theorem), the use of a monotonicity-preserving limiter reduces the convergence order of any algorithm to 1, so we do not include any test cases that use `T_lim!` and `lim!`. + +The test cases we use for this analysis are: + - `ark_analytic`, which uses a nonlinear `T_exp!` and a linear `T_imp!` + - `ark_analytic_sys` and `ark_onewaycouple_mri`, which use a linear `T_imp!` + - `ark_analytic_nonlin`, which uses a nonlinear `T_imp!` + - `1d_heat_equation` and `2d_heat_equation`, which use a nonlinear `T_exp!` and `dss!`, where the spatial discretization is implemented using `ClimaCore` ```@example include("report_gen.jl") ``` - -Plots for `IMEXAlgorithm`: - ![](output/convergence_ark_analytic_nonlin_imex_algorithms.png) - ![](output/convergence_ark_analytic_sys_imex_algorithms.png) ![](output/convergence_ark_analytic_imex_algorithms.png) + ![](output/convergence_ark_analytic_sys_imex_algorithms.png) ![](output/convergence_ark_onewaycouple_mri_imex_algorithms.png) + ![](output/convergence_ark_analytic_nonlin_imex_algorithms.png) ![](output/convergence_1d_heat_equation_imex_algorithms.png) ![](output/convergence_2d_heat_equation_imex_algorithms.png) +## Errors with a Limiter + +In order to verify the correctness of our algorithms with a limiter, we recreate Table 1 from ["Optimization-based limiters for the spectral element method" by Guba et al.](https://www.sciencedirect.com/science/article/pii/S0021999114001491) This involves running the `horizontal_deformational_flow` test case (from ["A standard test case suite for two-dimensional linear transport on the sphere" by Lauritzen et al.](https://gmd.copernicus.org/articles/5/887/2012/gmd-5-887-2012.pdf)) with and without a limiter, and also with and without hyperdiffusion. This test case uses a limited tendency `T_lim!` (which consists of advection and, optionally, hyperdiffusion), along with `dss!` and `lim!`. The spatial discretization is implemented using `ClimaCore`. Since this analysis is relatively expensive to run, we only check the results for `SSP333` and `ARS343`. Note that it is possible to limit undershoots and overshoots to 0 (up to floating-point roundoff error) when using the `SSPConstrained` `SSP333`, but not when using the `Unconstrained` `ARS343`. + +```@example +using ClimaTimeSteppers # hide +ENV["GKSwstype"] = "nul" # hide +include(joinpath(@__DIR__, "..", "plotting_utils.jl")) # hide +include(joinpath(pkgdir(ClimaTimeSteppers), "test", "problems.jl")) # hide + +# This also runs with num_steps = 1000, but with larger under/overshoots; 4800 +# is the value used in the paper. +limiter_summary(Float64, [SSP333(), ARS343()], horizontal_deformational_flow_test, 4800) +``` + +Plots of the tracer specific humidities that were used to compute this table are shown below. + ![](output/limiter_summary_SSP333.png) + ![](output/limiter_summary_ARS343.png) + ## References - [Example Programs for ARK ode (SUNDIALS)](http://runge.math.smu.edu/ARKode_example.pdf) + - ["Optimization-based limiters for the spectral element method" by Guba et al.](https://www.sciencedirect.com/science/article/pii/S0021999114001491) + - ["A standard test case suite for two-dimensional linear transport on the sphere" by Lauritzen et al.](https://gmd.copernicus.org/articles/5/887/2012/gmd-5-887-2012.pdf) diff --git a/docs/src/plotting_utils.jl b/docs/src/plotting_utils.jl index c16ce43c..1ebc7822 100644 --- a/docs/src/plotting_utils.jl +++ b/docs/src/plotting_utils.jl @@ -1,7 +1,9 @@ -import Plots +import Plots, Markdown +using ClimaCorePlots using Distributions: quantile, TDist using Printf: @sprintf using LaTeXStrings: latexstring +using PrettyTables: pretty_table, ft_printf """ predicted_convergence_order(algorithm_name, ode_function) @@ -82,7 +84,7 @@ function convergence_order(dts, errs, confidence) return order, order_uncertainty end -function test_imex_algorithms( +function verify_convergence( title, algorithm_names, test_case, @@ -95,6 +97,8 @@ function test_imex_algorithms( full_history_algorithm_name = nothing, average_function = array -> norm(array) / sqrt(length(array)), average_function_str = "RMS", + only_endpoints = false, + verbose = false, ) (; test_name, t_end, linear_implicit, analytic_sol) = test_case prob = test_case.split_prob @@ -106,8 +110,11 @@ function test_imex_algorithms( analytic_sol else ref_alg = IMEXAlgorithm(numerical_reference_algorithm_name, newtons_method) + ref_alg_str = string(nameof(typeof(numerical_reference_algorithm_name))) ref_dt = t_end / numerical_reference_num_steps - solve(deepcopy(prob), ref_alg; dt = ref_dt, save_everystep = true) + verbose && + @info "Generating numerical reference solution for $test_name with $ref_alg_str (dt = $ref_dt)..." + solve(deepcopy(prob), ref_alg; dt = ref_dt, save_everystep = !only_endpoints) end cur_avg_err(u, t) = average_function(abs.(u .- ref_sol(t))) @@ -127,6 +134,7 @@ function test_imex_algorithms( cur_avg_err_str = "\\textrm{current}\\_$net_avg_err_str" linestyles = (:solid, :dash, :dot, :dashdot, :dashdotdot) + marker_kwargs = (; markershape = :circle, markeralpha = 0.5, markerstrokewidth = 0) plot_kwargs = (; legendposition = :outerright, legendtitlefontpointsize = 8, @@ -170,16 +178,18 @@ function test_imex_algorithms( predicted_order = predicted_convergence_order(algorithm_name, prob.f) linestyle = linestyles[(predicted_order - 1) % length(linestyles) + 1] + verbose && @info "Running $test_name with $alg_str..." plot1_net_avg_errs = map(plot1_dts) do plot1_dt cur_avg_errs = solve( deepcopy(prob), alg; dt = plot1_dt, - save_everystep = true, + save_everystep = !only_endpoints, save_func = cur_avg_err, kwargshandle = DiffEqBase.KeywordArgSilent, ).u + verbose && @info "RMS_error(dt = $plot1_dt) = $(average_function(cur_avg_errs))" return average_function(cur_avg_errs) end order, order_uncertainty = convergence_order(plot1_dts, plot1_net_avg_errs, order_confidence_percent / 100) @@ -190,22 +200,14 @@ function test_imex_algorithms( else plot1_label = "$alg_str: \$$order_str\$" end + verbose && @info "Order = $order ± $order_uncertainty" if abs(order - predicted_order) > order_uncertainty @warn "Predicted order outside error bars for $alg_str ($test_name)" end if order_uncertainty > predicted_order / 10 @warn "Order uncertainty too large for $alg_str ($test_name)" end - Plots.plot!( - plot1, - plot1_dts, - plot1_net_avg_errs; - label = latexstring(plot1_label), - markershape = :circle, - markeralpha = 0.5, - markerstrokewidth = 0, - linestyle, - ) + Plots.plot!(plot1, plot1_dts, plot1_net_avg_errs; label = latexstring(plot1_label), linestyle, marker_kwargs...) # Remove all 0s from plot2_cur_avg_errs because they cannot be plotted on a # logarithmic scale. Record the extrema of plot2_cur_avg_errs to set ylim. @@ -213,7 +215,7 @@ function test_imex_algorithms( deepcopy(prob), alg; dt = default_dt, - save_everystep = true, + save_everystep = !only_endpoints, save_func = cur_avg_sol_and_err, kwargshandle = DiffEqBase.KeywordArgSilent, ) @@ -231,6 +233,7 @@ function test_imex_algorithms( plot2_cur_avg_sols; label = latexstring("$alg_str: \$$(si_str(plot2_net_avg_sol))\$"), linestyle, + (only_endpoints ? marker_kwargs : (;))..., ) Plots.plot!( plot2b, @@ -238,6 +241,7 @@ function test_imex_algorithms( plot2_cur_avg_errs; label = latexstring("$alg_str: \$$(si_str(plot2_net_avg_err))\$"), linestyle, + (only_endpoints ? marker_kwargs : (;))..., ) end @@ -252,7 +256,7 @@ function test_imex_algorithms( deepcopy(prob), history_alg; dt = default_dt, - save_everystep = true, + save_everystep = !only_endpoints, save_func = (u, t) -> u .- ref_sol(t), kwargshandle = DiffEqBase.KeywordArgSilent, ) @@ -285,7 +289,6 @@ function test_imex_algorithms( \$$net_avg_sol_str = $net_avg_sol_def_str\$ \$$net_avg_err_str = $net_avg_err_def_str\$""" if !isnothing(numerical_reference_algorithm_name) - ref_alg_name = string(nameof(typeof(numerical_reference_algorithm_name))) ref_cur_avg_errs = map(ref_sol.u, ref_sol.t) do u, t average_function(abs.(u .- analytic_sol(t))) end @@ -294,7 +297,7 @@ function test_imex_algorithms( avg_def_str("|Y_{ref}(t)[\\textrm{index}] - Y_{analytic}(t)[\\textrm{index}]|", "\\textrm{index}") ref_net_avg_err_def_str = avg_def_str(ref_cur_avg_err_def_str, "t") footnote = "$footnote\n\n\nNote: The \"reference solution\" \$Y_{ref}\$ was \ - computed using $ref_alg_name with\n\$dt = $(pow_str(ref_dt)),\\ \ + computed using $ref_alg_str with\n\$dt = $(pow_str(ref_dt)),\\ \ \\textrm{and}\\ $ref_net_avg_err_def_str = \ $(si_str(ref_net_avg_err))\$" end @@ -321,3 +324,95 @@ function test_imex_algorithms( file_suffix = lowercase(replace(test_name * ' ' * title, " " => "_")) Plots.savefig(plot, joinpath("output", "convergence_$file_suffix.png")) end + +# Generates the Table 1 from +# "Optimization-based limiters for the spectral element method" by Guba et al., +# and also plots the values used to generate the table. +function limiter_summary(::Type{FT}, algorithm_names, test_case_type, num_steps) where {FT} + to_title(name) = titlecase(replace(string(name), '_' => ' ')) + table_rows = [] + for algorithm_name in algorithm_names + alg_str = string(nameof(typeof(algorithm_name))) + plots = [] + plot_kwargs = (; + clims = (0, 1), + color = :diverging_rainbow_bgymr_45_85_c67_n256, + colorbar = false, + guide = "", + margin = 10Plots.px, + ) + for use_limiter in (false, true), use_hyperdiffusion in (false, true) + test_case = test_case_type(FT; use_limiter, use_hyperdiffusion) + prob = test_case.split_prob + dt = test_case.t_end / num_steps + algorithm = IMEXAlgorithm(algorithm_name, NewtonsMethod()) + solution = solve(deepcopy(prob), algorithm; dt).u + initial_q = solution[1].ρq ./ solution[1].ρ + final_q = solution[end].ρq ./ solution[end].ρ + names = propertynames(initial_q) + + if isempty(plots) + for name in names + push!(plots, Plots.plot(initial_q.:($name); plot_kwargs..., title = to_title(name))) + end + end + for name in names + push!(plots, Plots.plot(final_q.:($name); plot_kwargs..., title = "")) + end + + for name in names + ϕ₀ = initial_q.:($name) + ϕ = final_q.:($name) + Δϕ₀ = maximum(ϕ₀) - minimum(ϕ₀) + ϕ_error = ϕ .- ϕ₀ + table_row = [ + alg_str;; + string(use_limiter);; + string(use_hyperdiffusion);; + to_title(name);; + max(0, -(minimum(ϕ) - minimum(ϕ₀)) / Δϕ₀);; + max(0, (maximum(ϕ) - maximum(ϕ₀)) / Δϕ₀);; + map(p -> norm(ϕ_error, p) / norm(ϕ₀, p), (1, 2, Inf))... + ] + push!(table_rows, table_row) + end + end + colorbar_plot = Plots.scatter( + [0]; + plot_kwargs..., + colorbar = true, + framestyle = :none, + legend_position = :none, + margin = 0Plots.px, + markeralpha = 0, + zcolor = [0], + ) + plot = Plots.plot( + plots..., + colorbar_plot; + layout = (Plots.@layout [Plots.grid(5, 3) a{0.1w}]), + plot_title = "Tracer specific humidity for $alg_str (Initial, \ + Final, Final w/ Hyperdiffusion, Final w/ Limiter, \ + Final w/ Hyperdiffusion & Limiter)", + size = (1600, 2000), + ) + Plots.savefig(plot, joinpath("output", "limiter_summary_$alg_str.png")) + end + table = pretty_table( + vcat(table_rows...); + header = [ + "Algorithm", + "Limiter", + "Hyperdiffusion", + "Tracer Name", + "Max Undershoot", + "Max Overshoot", + "1-Norm Error", + "2-Norm Error", + "∞-Norm Error", + ], + body_hlines = collect(3:3:(length(table_rows) - 1)), + formatters = ft_printf("%.4e"), + ) + println(table) +end diff --git a/src/solvers/imex_ssp.jl b/src/solvers/imex_ssp.jl index f8004c62..08fd51c1 100644 --- a/src/solvers/imex_ssp.jl +++ b/src/solvers/imex_ssp.jl @@ -80,7 +80,7 @@ function step_u!(integrator, cache::IMEXSSPRKCache) elseif !iszero(β[i - 1]) if !isnothing(T_lim!) @. U_lim = U_exp + dt * T_lim - lim!(U_lim, p, U_exp) + lim!(U_lim, p, t_exp, U_exp) @. U_exp = U_lim end if !isnothing(T_exp!) @@ -151,7 +151,7 @@ function step_u!(integrator, cache::IMEXSSPRKCache) if !iszero(β[s]) if !isnothing(T_lim!) @. U_lim = U_exp + dt * T_lim - lim!(U_lim, p, U_exp) + lim!(U_lim, p, t_final, U_exp) @. U_exp = U_lim end if !isnothing(T_exp!) diff --git a/test/problems.jl b/test/problems.jl index a93dbc85..5be746b5 100644 --- a/test/problems.jl +++ b/test/problems.jl @@ -534,6 +534,308 @@ function climacore_1Dheat_test_cts(::Type{FT}) where {FT} ) end +# Monkey patch for ClimaCore so that ρq .* current_wind_vector works +ClimaCore.RecursiveApply.rmap(fn::F, X::Tuple, Y) where {F} = + ClimaCore.RecursiveApply.tuplemap(x -> ClimaCore.RecursiveApply.rmap(fn, x, Y), X) +ClimaCore.RecursiveApply.rmap(fn::F, X::NamedTuple{names}, Y) where {F, names} = + NamedTuple{names}(ClimaCore.RecursiveApply.rmap(fn, Tuple(X), Y)) + +# "Dynamical Core Model Intercomparison Project (DCMIP) Test Case Document" by +# Ullrich et al., Section 1.1 +# (http://www-personal.umich.edu/~cjablono/DCMIP-2012_TestCaseDocument_v1.7.pdf) +# Implemented in flux form, with an optional limiter and hyperdiffusion. +# TODO: Use this as an integration test. +function deformational_flow_test(::Type{FT}; use_limiter = true, use_hyperdiffusion = true) where {FT} + # Table III + # Note: the paper uses "a" in place of "R" + R = FT(6371220) # radius of Earth [m] + g = FT(9.80616) # gravitational acceleration [m/s^2] + R_d = FT(287.0) # gas constant for dry air [J/kg/K] + + # Table IX + # Note: the paper specifies that λ_c1 = 5π/6 ∈ [0, 2π), and that + # λ_c2 = 7π/6 ∈ [0, 2π), whereas we use λ_c1, λ_c2 ∈ [-π, π). + z_top = FT(12000) # altitude at model top [m] + p_top = FT(25494.4) # pressure at model top [Pa] + T_0 = FT(300) # isothermal atmospheric temperature [K] + p_0 = FT(100000) # reference pressure [Pa] + τ = 60 * 60 * 24 * FT(12) # period of motion (12 days) [s] + ω_0 = 23000 * FT(π) / τ # maximum vertical pressure velocity [Pa/s] + b = FT(0.2) # normalized pressure depth of divergent layer + λ_c1 = -FT(π) / 6 # initial longitude of first tracer + λ_c2 = FT(π) / 6 # initial longitude of second tracer + φ_c = FT(0) # initial latitude of tracers + z_c = FT(5000) # initial altitude of tracers [m] + R_t = R / 2 # horizontal half-width of tracers [m] + Z_t = FT(1000) # vertical half-width of tracers [m] + + # scale height [m] (Equation 3) + H = R_d * T_0 / g + + # hyperviscosity coefficient [m^4] (specified in the limiter paper) + D₄ = FT(6.6e14) + + centers = ClimaCore.Geometry.LatLongZPoint.(rad2deg(φ_c), rad2deg.((λ_c1, λ_c2)), FT(0)) + + # custom discretization (paper's discretization results in a very slow test) + vert_nelems = 10 + horz_nelems = 4 + horz_npoly = 3 + + vert_domain = ClimaCore.Domains.IntervalDomain( + ClimaCore.Geometry.ZPoint(FT(0)), + ClimaCore.Geometry.ZPoint(z_top); + boundary_names = (:bottom, :top), + ) + vert_mesh = ClimaCore.Meshes.IntervalMesh(vert_domain, nelems = vert_nelems) + vert_cent_space = ClimaCore.Spaces.CenterFiniteDifferenceSpace(vert_mesh) + + horz_domain = ClimaCore.Domains.SphereDomain(R) + horz_mesh = ClimaCore.Meshes.EquiangularCubedSphere(horz_domain, horz_nelems) + horz_topology = ClimaCore.Topologies.Topology2D(horz_mesh) + horz_quad = ClimaCore.Spaces.Quadratures.GLL{horz_npoly + 1}() + horz_space = ClimaCore.Spaces.SpectralElementSpace2D(horz_topology, horz_quad) + + cent_space = ClimaCore.Spaces.ExtrudedFiniteDifferenceSpace(horz_space, vert_cent_space) + cent_coords = ClimaCore.Fields.coordinate_field(cent_space) + face_space = ClimaCore.Spaces.FaceExtrudedFiniteDifferenceSpace(cent_space) + + # initial density (Equation 8) + cent_ρ = @. p_0 / (R_d * T_0) * exp(-cent_coords.z / H) + + # initial tracer concentrations (Equations 28--35) + cent_q = map(cent_coords) do coord + z = coord.z + φ = deg2rad(coord.lat) + + ds = map(centers) do center + r = ClimaCore.Geometry.great_circle_distance(coord, center, horz_space.global_geometry) + return min(1, (r / R_t)^2 + ((z - z_c) / Z_t)^2) + end + in_slot = z > z_c && φ_c - FT(0.125) < φ < φ_c + FT(0.125) + + q1 = (1 + cos(FT(π) * ds[1])) / 2 + (1 + cos(FT(π) * ds[2])) / 2 + q2 = FT(0.9) - FT(0.8) * q1^2 + q3 = (ds[1] < FT(0.5) || ds[2] < FT(0.5)) && !in_slot ? FT(1) : FT(0.1) + q4 = 1 - FT(0.3) * (q1 + q2 + q3) + q5 = FT(1) + return (; q1, q2, q3, q4, q5) + end + + init_state = ClimaCore.Fields.FieldVector(; cent_ρ, cent_ρq = cent_ρ .* cent_q) + + # current wind vector (Equations 15--26) + current_cent_wind_vector = ClimaCore.Fields.Field(ClimaCore.Geometry.UVWVector{FT}, cent_space) + current_face_wind_vector = ClimaCore.Fields.Field(ClimaCore.Geometry.UVWVector{FT}, face_space) + function wind_vector(coord, ρ, t) + z = coord.z + φ = deg2rad(coord.lat) + λ = deg2rad(coord.long) + + p = p_0 * exp(-g * z / (R_d * T_0)) # initial pressure (Equation 1) + λ′ = λ - 2 * FT(π) * t / τ + k = 10 * R / τ + + u_a = k * sin(λ′)^2 * sin(2 * φ) * cos(FT(π) * t / τ) + 2 * FT(π) * R / τ * cos(φ) + u_d = + ω_0 * R / (b * p_top) * + cos(λ′) * + cos(φ)^2 * + cos(2 * FT(π) * t / τ) * + (-exp((p - p_0) / (b * p_top)) + exp((p_top - p) / (b * p_top))) + u = u_a + u_d + v = k * sin(2 * λ′) * cos(φ) * cos(FT(π) * t / τ) + s = 1 + exp((p_top - p_0) / (b * p_top)) - exp((p - p_0) / (b * p_top)) - exp((p_top - p) / (b * p_top)) + ω = ω_0 * sin(λ′) * cos(φ) * cos(2 * FT(π) * t / τ) * s + w = -ω / (g * ρ) + + return ClimaCore.Geometry.UVWVector(u, v, w) + end + + horz_div = ClimaCore.Operators.Divergence() + horz_wdiv = ClimaCore.Operators.WeakDivergence() + horz_grad = ClimaCore.Operators.Gradient() + cent_χ = similar(cent_q) + function T_lim!(tendency, state, _, t) + @. current_cent_wind_vector = wind_vector(cent_coords, state.cent_ρ, t) + @. tendency.cent_ρ = -horz_div(state.cent_ρ * current_cent_wind_vector) + @. tendency.cent_ρq = -horz_div(state.cent_ρq * current_cent_wind_vector) + use_hyperdiffusion || return nothing + @. cent_χ = horz_wdiv(horz_grad(state.cent_ρq / state.cent_ρ)) + ClimaCore.Spaces.weighted_dss!(cent_χ) + @. tendency.cent_ρq += -D₄ * horz_wdiv(state.cent_ρ * horz_grad(cent_χ)) + return nothing + end + + limiter = ClimaCore.Limiters.QuasiMonotoneLimiter(cent_q; rtol = FT(0)) + function lim!(state, _, t, ref_state) + use_limiter || return nothing + ClimaCore.Limiters.compute_bounds!(limiter, ref_state.cent_ρq, ref_state.cent_ρ) + ClimaCore.Limiters.apply_limiter!(state.cent_ρq, state.cent_ρ, limiter) + return nothing + end + + vert_div = ClimaCore.Operators.DivergenceF2C() + vert_interp = ClimaCore.Operators.InterpolateC2F( + top = ClimaCore.Operators.Extrapolate(), + bottom = ClimaCore.Operators.Extrapolate(), + ) + function T_exp!(tendency, state, _, t) + @. current_face_wind_vector = wind_vector(face_coords, vert_interp(state.cent_ρ), t) + @. tendency.cent_ρ = -vert_div(vert_interp(state.cent_ρ) * current_face_wind_vector) + @. tendency.cent_ρq = -vert_div(vert_interp(state.cent_ρq) * current_face_wind_vector) + end + + function dss!(state, _, t) + ClimaCore.Spaces.weighted_dss!(state.q) + end + + function analytic_sol(t) + t ∈ (0, τ) || error("Analytic solution only defined at start and end") + return copy(init_state) + end + + tendency_func = ClimaODEFunction(; T_lim!, T_exp!, lim!, dss!) + split_tendency_func = tendency_func + make_prob(func) = ODEProblem(func, init_state, (FT(0), τ), nothing) + IntegratorTestCase( + "Deformational Flow", + false, + τ, + analytic_sol, + make_prob(tendency_func), + make_prob(split_tendency_func), + ) +end + +# "A standard test case suite for two-dimensional linear transport on the +# sphere" by Lauritzen et al. +# (https://gmd.copernicus.org/articles/5/887/2012/gmd-5-887-2012.pdf) +# Implemented in flux form, with an optional limiter and hyperdiffusion. +function horizontal_deformational_flow_test(::Type{FT}; use_limiter = true, use_hyperdiffusion = true) where {FT} + # constants (using the same notation as deformational_flow_test) + R = FT(6371220) # radius of Earth [m] + τ = 60 * 60 * 24 * FT(12) # period of motion (12 days) [s] + λ_c1 = -FT(π) / 6 # initial longitude of first tracer + λ_c2 = FT(π) / 6 # initial longitude of second tracer + φ_c = FT(0) # initial latitude of tracers + R_t = R / 2 # horizontal half-width of tracers [m] + D₄ = FT(6.6e14) # hyperviscosity coefficient [m^4] (specified in the limiter paper) + + centers = ClimaCore.Geometry.LatLongPoint.(rad2deg(φ_c), rad2deg.((λ_c1, λ_c2))) + + # 1.5° resolution on the equator: 360° / (4 * nelems * npoly) = 1.5° + nelems = 20 + npoly = 3 + + domain = ClimaCore.Domains.SphereDomain(R) + mesh = ClimaCore.Meshes.EquiangularCubedSphere(domain, nelems) + topology = ClimaCore.Topologies.Topology2D(mesh) + quad = ClimaCore.Spaces.Quadratures.GLL{npoly + 1}() + space = ClimaCore.Spaces.SpectralElementSpace2D(topology, quad) + coords = ClimaCore.Fields.coordinate_field(space) + + # initial conditions (Section 2.2) + ρ = ones(space) + q = map(coords) do coord + φ = deg2rad(coord.lat) + λ = deg2rad(coord.long) + + hs = map(centers) do center + center′ = ClimaCore.Geometry.CartesianPoint(center, space.global_geometry) + coord′ = ClimaCore.Geometry.CartesianPoint(coord, space.global_geometry) + dist_squared = (coord′.x1 - center′.x1)^2 + (coord′.x2 - center′.x2)^2 + (coord′.x3 - center′.x3)^2 + # Note: the paper doesn't divide by R^2, which only works if R = 1 + return FT(0.95) * exp(-5 * dist_squared / R^2) + end + gaussian_hills = hs[1] + hs[2] + rs = map(centers) do center + return ClimaCore.Geometry.great_circle_distance(coord, center, space.global_geometry) + end + cosine_bells = if rs[1] < R_t + FT(0.1) + FT(0.9) * (1 + cos(FT(π) * rs[1] / R_t)) / 2 + elseif rs[2] < R_t + FT(0.1) + FT(0.9) * (1 + cos(FT(π) * rs[2] / R_t)) / 2 + else + FT(0.1) + end + slotted_cylinders = + if ( + (rs[1] <= R_t && abs(λ - λ_c1) >= R_t / 6R) || + (rs[2] <= R_t && abs(λ - λ_c2) >= R_t / 6R) || + (rs[1] <= R_t && abs(λ - λ_c1) < R_t / 6R && φ - φ_c < -5R_t / 12R) || + (rs[2] <= R_t && abs(λ - λ_c2) < R_t / 6R && φ - φ_c > 5R_t / 12R) + ) + FT(1) + else + FT(0.1) + end + + return (; gaussian_hills, cosine_bells, slotted_cylinders) + end + init_state = ClimaCore.Fields.FieldVector(; ρ, ρq = ρ .* q) + + # current wind vector (Section 2.3) + current_wind_vector = ClimaCore.Fields.Field(ClimaCore.Geometry.UVVector{FT}, space) + function wind_vector(coord, t) + φ = deg2rad(coord.lat) + λ = deg2rad(coord.long) + + λ′ = λ - 2 * FT(π) * t / τ + k = 10 * R / τ + + u = k * sin(λ′)^2 * sin(2 * φ) * cos(FT(π) * t / τ) + 2 * FT(π) * R / τ * cos(φ) + v = k * sin(2 * λ′) * cos(φ) * cos(FT(π) * t / τ) + + return ClimaCore.Geometry.UVVector(u, v) + end + + div = ClimaCore.Operators.Divergence() + wdiv = ClimaCore.Operators.WeakDivergence() + grad = ClimaCore.Operators.Gradient() + χ = similar(q) + function T_lim!(tendency, state, _, t) + @. current_wind_vector = wind_vector(coords, t) + @. tendency.ρ = -div(state.ρ * current_wind_vector) + @. tendency.ρq = -div(state.ρq * current_wind_vector) + use_hyperdiffusion || return nothing + @. χ = wdiv(grad(state.ρq / state.ρ)) + ClimaCore.Spaces.weighted_dss!(χ) + @. tendency.ρq += -D₄ * wdiv(state.ρ * grad(χ)) + return nothing + end + + limiter = ClimaCore.Limiters.QuasiMonotoneLimiter(q; rtol = FT(0)) + function lim!(state, _, t, ref_state) + use_limiter || return nothing + ClimaCore.Limiters.compute_bounds!(limiter, ref_state.ρq, ref_state.ρ) + ClimaCore.Limiters.apply_limiter!(state.ρq, state.ρ, limiter) + return nothing + end + + function dss!(state, _, t) + ClimaCore.Spaces.weighted_dss!(state.ρ) + ClimaCore.Spaces.weighted_dss!(state.ρq) + end + + function analytic_sol(t) + t ∈ (0, τ) || error("Analytic solution only defined at start and end") + return copy(init_state) + end + + tendency_func = ClimaODEFunction(; T_lim!, lim!, dss!) + split_tendency_func = tendency_func + make_prob(func) = ODEProblem(func, init_state, (FT(0), τ), nothing) + IntegratorTestCase( + "Horizontal Deformational Flow", + false, + τ, + analytic_sol, + make_prob(tendency_func), + make_prob(split_tendency_func), + ) +end + function all_test_cases(::Type{FT}) where {FT} return [ ark_analytic_nonlin_test_cts(FT),