perf: added matrix operations and multithreading
Using matrix operations wherever possible Added support for multithreading Moved simulation loop into BTCS to minimize memory allocation Switched to Tridiagonal Coefficient Matrix [skip cli]
This commit is contained in:
parent
a064f9de24
commit
3c080c7149
@ -133,7 +133,7 @@ def main(tolerance, runs, silent, no_clean, precompile):
|
|||||||
|
|
||||||
# Print results
|
# Print results
|
||||||
if not silent: print("\n----- Benchmark Results -----")
|
if not silent: print("\n----- Benchmark Results -----")
|
||||||
print(f"Parameters: Tolerance = {tolerance}, Runs = {runs}")
|
print(f"Parameters: Tolerance = {tolerance}, Runs = {runs}, Precompile = {precompile}")
|
||||||
for name in sorted(results_dict):
|
for name in sorted(results_dict):
|
||||||
print(results_dict[name])
|
print(results_dict[name])
|
||||||
if pass_all:
|
if pass_all:
|
||||||
|
|||||||
@ -5,6 +5,8 @@
|
|||||||
|
|
||||||
using LinearAlgebra
|
using LinearAlgebra
|
||||||
using SparseArrays
|
using SparseArrays
|
||||||
|
using Base.Threads
|
||||||
|
using CUDA
|
||||||
|
|
||||||
include("../Boundary.jl")
|
include("../Boundary.jl")
|
||||||
include("../Grid.jl")
|
include("../Grid.jl")
|
||||||
@ -13,6 +15,10 @@ function calcAlphaIntercell(alpha1::T, alpha2::T) where {T}
|
|||||||
return 2 / ((1 / alpha1) + (1 / alpha2))
|
return 2 / ((1 / alpha1) + (1 / alpha2))
|
||||||
end
|
end
|
||||||
|
|
||||||
|
function calcAlphaIntercell(alpha1::Matrix{T}, alpha2::Matrix{T}) where {T}
|
||||||
|
return 2 ./ ((1 ./ alpha1) .+ (1 ./ alpha2))
|
||||||
|
end
|
||||||
|
|
||||||
function calcBoundaryCoeffConstant(alpha_center::T, alpha_side::T, sx::T) where {T}
|
function calcBoundaryCoeffConstant(alpha_center::T, alpha_side::T, sx::T) where {T}
|
||||||
alpha = calcAlphaIntercell(alpha_center, alpha_side)
|
alpha = calcAlphaIntercell(alpha_center, alpha_side)
|
||||||
centerCoeff = 1 + sx * (alpha + 2 * alpha_center)
|
centerCoeff = 1 + sx * (alpha + 2 * alpha_center)
|
||||||
@ -28,47 +34,37 @@ function calcBoundaryCoeffClosed(alpha_center::T, alpha_side::T, sx::T) where {T
|
|||||||
end
|
end
|
||||||
|
|
||||||
# creates coefficient matrix for next time step from alphas in x-direction
|
# creates coefficient matrix for next time step from alphas in x-direction
|
||||||
function createCoeffMatrix(alpha::Matrix{T}, bcLeft::Vector{BoundaryElement{T}}, bcRight::Vector{BoundaryElement{T}}, numCols::Int, rowIndex::Int, sx::T) where {T}
|
function createCoeffMatrix(alpha::Matrix{T}, alpha_left::Vector{T}, alpha_right::Vector{T}, bcLeft::Vector{BoundaryElement{T}}, bcRight::Vector{BoundaryElement{T}}, numCols::Int, rowIndex::Int, sx::T)::Tridiagonal{T} where {T}
|
||||||
numCols = max(numCols, 2)
|
# Precompute boundary condition type check for efficiency
|
||||||
cm = spzeros(T, numCols, numCols)
|
bcTypeLeft = getType(bcLeft[rowIndex])
|
||||||
|
|
||||||
# left column
|
# Determine left side boundary coefficients based on boundary condition
|
||||||
if getType(bcLeft[rowIndex]) == CONSTANT
|
centerCoeffTop, rightCoeffTop = if bcTypeLeft == CONSTANT
|
||||||
centerCoeffTop, rightCoeffTop = calcBoundaryCoeffConstant(alpha[rowIndex, 1], alpha[rowIndex, 2], sx)
|
calcBoundaryCoeffConstant(alpha[rowIndex, 1], alpha[rowIndex, 2], sx)
|
||||||
cm[1, 1] = centerCoeffTop
|
elseif bcTypeLeft == CLOSED
|
||||||
cm[1, 2] = rightCoeffTop
|
calcBoundaryCoeffClosed(alpha[rowIndex, 1], alpha[rowIndex, 2], sx)
|
||||||
elseif getType(bcLeft[rowIndex]) == CLOSED
|
|
||||||
centerCoeffTop, rightCoeffTop = calcBoundaryCoeffClosed(alpha[rowIndex, 1], alpha[rowIndex, 2], sx)
|
|
||||||
cm[1, 1] = centerCoeffTop
|
|
||||||
cm[1, 2] = rightCoeffTop
|
|
||||||
else
|
else
|
||||||
error("Undefined Boundary Condition Type somewhere on Left or Top!")
|
error("Undefined Boundary Condition Type on Left!")
|
||||||
end
|
end
|
||||||
|
|
||||||
# inner columns
|
# Precompute boundary condition type check for efficiency
|
||||||
@inbounds for i in 2:(numCols-1)
|
bcTypeRight = getType(bcRight[rowIndex])
|
||||||
alpha_left_here = calcAlphaIntercell(alpha[rowIndex, i-1], alpha[rowIndex, i])
|
|
||||||
alpha_here_right = alpha[rowIndex, i-1] == alpha[rowIndex, i+1] ? alpha_left_here : calcAlphaIntercell(alpha[rowIndex, i], alpha[rowIndex, i+1]) # calcAlphaIntercell is symmetric, so we can use it for both directions
|
|
||||||
|
|
||||||
cm[i, i-1] = -sx * alpha_left_here
|
# Determine right side boundary coefficients based on boundary condition
|
||||||
cm[i, i] = 1 + sx * (alpha_here_right + alpha_left_here)
|
centerCoeffBottom, leftCoeffBottom = if bcTypeRight == CONSTANT
|
||||||
cm[i, i+1] = -sx * alpha_here_right
|
calcBoundaryCoeffConstant(alpha[rowIndex, numCols], alpha[rowIndex, numCols-1], sx)
|
||||||
end
|
elseif bcTypeRight == CLOSED
|
||||||
|
calcBoundaryCoeffClosed(alpha[rowIndex, numCols], alpha[rowIndex, numCols-1], sx)
|
||||||
# right column
|
|
||||||
if getType(bcRight[rowIndex]) == CONSTANT
|
|
||||||
centerCoeffBottom, leftCoeffBottom = calcBoundaryCoeffConstant(alpha[rowIndex, numCols], alpha[rowIndex, numCols-1], sx)
|
|
||||||
cm[numCols, numCols-1] = leftCoeffBottom
|
|
||||||
cm[numCols, numCols] = centerCoeffBottom
|
|
||||||
elseif getType(bcRight[rowIndex]) == CLOSED
|
|
||||||
centerCoeffBottom, leftCoeffBottom = calcBoundaryCoeffClosed(alpha[rowIndex, numCols], alpha[rowIndex, numCols-1], sx)
|
|
||||||
cm[numCols, numCols-1] = leftCoeffBottom
|
|
||||||
cm[numCols, numCols] = centerCoeffBottom
|
|
||||||
else
|
else
|
||||||
error("Undefined Boundary Condition Type somewhere on Right or Bottom!")
|
error("Undefined Boundary Condition Type on Right!")
|
||||||
end
|
end
|
||||||
|
|
||||||
return cm
|
dl = [-sx .* alpha_left; leftCoeffBottom]
|
||||||
|
du = [rightCoeffTop; -sx .* alpha_right]
|
||||||
|
d = [centerCoeffTop; 1 .+ sx .* (alpha_right + alpha_left); centerCoeffBottom]
|
||||||
|
alpha_diagonal = Tridiagonal(dl, d, du)
|
||||||
|
|
||||||
|
return alpha_diagonal
|
||||||
end
|
end
|
||||||
|
|
||||||
|
|
||||||
@ -86,9 +82,9 @@ function calcExplicitConcentrationsBoundaryConstant(conc_center::T, conc_bc::T,
|
|||||||
sy * alpha_center * conc_bc
|
sy * alpha_center * conc_bc
|
||||||
end
|
end
|
||||||
|
|
||||||
function createSolutionVector(concentrations::Matrix{T}, alphaX::Matrix{T}, alphaY::Matrix{T}, bcLeft::Vector{BoundaryElement{T}}, bcRight::Vector{BoundaryElement{T}}, bcTop::Vector{BoundaryElement{T}}, bcBottom::Vector{BoundaryElement{T}}, length::Int, rowIndex::Int, sx::T, sy::T) where {T}
|
function writeSolutionVector!(sv::Vector{T}, concentrations::Matrix{T}, alphaX::Matrix{T}, alphaY::Matrix{T}, bcLeft::Vector{BoundaryElement{T}}, bcRight::Vector{BoundaryElement{T}}, bcTop::Vector{BoundaryElement{T}}, bcBottom::Vector{BoundaryElement{T}}, rowIndex::Int, sx::T, sy::T) where {T}
|
||||||
numRows = size(concentrations, 1)
|
numRows = size(concentrations, 1)
|
||||||
sv = Vector{T}(undef, length)
|
length = size(sv, 1)
|
||||||
|
|
||||||
# Inner rows
|
# Inner rows
|
||||||
if rowIndex > 1 && rowIndex < numRows
|
if rowIndex > 1 && rowIndex < numRows
|
||||||
@ -136,33 +132,28 @@ function createSolutionVector(concentrations::Matrix{T}, alphaX::Matrix{T}, alph
|
|||||||
if getType(bcRight[rowIndex]) == CONSTANT
|
if getType(bcRight[rowIndex]) == CONSTANT
|
||||||
sv[end] += 2 * sx * alphaX[rowIndex, end] * getValue(bcRight[rowIndex])
|
sv[end] += 2 * sx * alphaX[rowIndex, end] * getValue(bcRight[rowIndex])
|
||||||
end
|
end
|
||||||
|
|
||||||
return sv
|
|
||||||
end
|
|
||||||
|
|
||||||
# solver for linear equation system; A corresponds to coefficient matrix, b to the solution vector
|
|
||||||
function LinearAlgebraAlgorithm(A::SparseMatrixCSC{T}, b::Vector{T}) where {T}
|
|
||||||
return A \ b
|
|
||||||
end
|
end
|
||||||
|
|
||||||
# BTCS solution for 1D grid
|
# BTCS solution for 1D grid
|
||||||
function BTCS_1D(grid::Grid{T}, bc::Boundary{T}, timestep::T, solverFunc::Function) where {T}
|
function BTCS_1D(grid::Grid{T}, bc::Boundary{T}, timestep::T) where {T}
|
||||||
length = grid.cols
|
length = getCols(grid)
|
||||||
sx = timestep / (grid.deltaCol * grid.deltaCol)
|
sx = timestep / (grid.deltaCol * grid.deltaCol)
|
||||||
|
|
||||||
b = Vector{T}(undef, length)
|
|
||||||
|
|
||||||
alpha = getAlphaX(grid)
|
alpha = getAlphaX(grid)
|
||||||
|
|
||||||
bcLeft = getBoundarySide(bc, LEFT)
|
bcLeft = getBoundarySide(bc, LEFT)
|
||||||
bcRight = getBoundarySide(bc, RIGHT)
|
bcRight = getBoundarySide(bc, RIGHT)
|
||||||
|
|
||||||
concentrations = grid.concentrations[]
|
concentrations::Matrix{T} = grid.concentrations[]
|
||||||
rowIndex = 1
|
rowIndex = 1
|
||||||
A = createCoeffMatrix(alpha, bcLeft, bcRight, length, rowIndex, sx)
|
|
||||||
@inbounds for i in 1:length
|
numCols = max(length, 2)
|
||||||
b[i] = concentrations[1, i]
|
|
||||||
end
|
alpha_left = calcAlphaIntercell(alpha[:, 1:(numCols-2)], alpha[:, 2:(numCols-1)])
|
||||||
|
alpha_right = calcAlphaIntercell(alpha[:, 2:(numCols-1)], alpha[:, 3:numCols])
|
||||||
|
|
||||||
|
A::Tridiagonal{T} = createCoeffMatrix(alpha, alpha_left[rowIndex, :], alpha_right[rowIndex, :], bcLeft, bcRight, length, rowIndex, sx)
|
||||||
|
b = concentrations[1, :]
|
||||||
|
|
||||||
if getType(getBoundarySide(bc, LEFT)[1]) == CONSTANT
|
if getType(getBoundarySide(bc, LEFT)[1]) == CONSTANT
|
||||||
b[1] += 2 * sx * alpha[1, 1] * bcLeft[1].value
|
b[1] += 2 * sx * alpha[1, 1] * bcLeft[1].value
|
||||||
@ -171,70 +162,90 @@ function BTCS_1D(grid::Grid{T}, bc::Boundary{T}, timestep::T, solverFunc::Functi
|
|||||||
b[length] += 2 * sx * alpha[1, length] * bcRight[1].value
|
b[length] += 2 * sx * alpha[1, length] * bcRight[1].value
|
||||||
end
|
end
|
||||||
|
|
||||||
concentrations_t1 = solverFunc(A, b)
|
concentrations_t1 = A \ b
|
||||||
|
|
||||||
@inbounds for j in 1:length
|
concentrations[1, :] = concentrations_t1
|
||||||
concentrations[1, j] = concentrations_t1[j]
|
|
||||||
end
|
|
||||||
|
|
||||||
setConcentrations!(grid, concentrations)
|
setConcentrations!(grid, concentrations)
|
||||||
end
|
end
|
||||||
|
|
||||||
# BTCS solution for 2D grid
|
# BTCS solution for 2D grid
|
||||||
function BTCS_2D(grid::Grid{T}, bc::Boundary{T}, timestep::T, solverFunc::Function, numThreads::Int) where {T}
|
function BTCS_2D(grid::Grid{T}, bc::Boundary{T}, alphaX_left::Matrix{T}, alphaX_right::Matrix{T}, alphaY_t_left::Matrix{T}, alphaY_t_right::Matrix{T}, timestep::T) where {T}
|
||||||
rowMax = grid.rows
|
rows = getRows(grid)
|
||||||
colMax = grid.cols
|
cols = getCols(grid)
|
||||||
sx = timestep / (2 * grid.deltaCol * grid.deltaCol)
|
sx = timestep / (2 * grid.deltaCol * grid.deltaCol)
|
||||||
sy = timestep / (2 * grid.deltaRow * grid.deltaRow)
|
sy = timestep / (2 * grid.deltaRow * grid.deltaRow)
|
||||||
|
|
||||||
concentrations_t1 = zeros(T, rowMax, colMax)
|
|
||||||
row_t1 = Vector{T}(undef, colMax)
|
|
||||||
|
|
||||||
alphaX = getAlphaX(grid)
|
alphaX = getAlphaX(grid)
|
||||||
alphaY = getAlphaY(grid)
|
alphaY = getAlphaY(grid)
|
||||||
|
|
||||||
|
alphaX_t = getAlphaX_t(grid)
|
||||||
|
alphaY_t = getAlphaY_t(grid)
|
||||||
|
|
||||||
|
concentrations = getConcentrations(grid)
|
||||||
|
concentrations_intermediate = similar(concentrations)
|
||||||
|
concentrations_t_task = Threads.@spawn copy(transpose(concentrations))
|
||||||
|
|
||||||
bcLeft = getBoundarySide(bc, LEFT)
|
bcLeft = getBoundarySide(bc, LEFT)
|
||||||
bcRight = getBoundarySide(bc, RIGHT)
|
bcRight = getBoundarySide(bc, RIGHT)
|
||||||
bcTop = getBoundarySide(bc, TOP)
|
bcTop = getBoundarySide(bc, TOP)
|
||||||
bcBottom = getBoundarySide(bc, BOTTOM)
|
bcBottom = getBoundarySide(bc, BOTTOM)
|
||||||
|
|
||||||
concentrations = grid.concentrations[]
|
localBs = [zeros(T, cols) for _ in 1:Threads.nthreads()]
|
||||||
|
Threads.@threads for i = 1:rows
|
||||||
|
localB = localBs[Threads.threadid()]
|
||||||
|
A::Tridiagonal{T} = createCoeffMatrix(alphaX, alphaX_left[i, :], alphaX_right[i, :], bcLeft, bcRight, cols, i, sx)
|
||||||
|
writeSolutionVector!(localB, concentrations, alphaX, alphaY, bcLeft, bcRight, bcTop, bcBottom, i, sx, sy)
|
||||||
|
|
||||||
@inbounds for i = 1:rowMax
|
concentrations_intermediate[i, :] = A \ localB
|
||||||
A = createCoeffMatrix(alphaX, bcLeft, bcRight, colMax, i, sx)
|
|
||||||
b = createSolutionVector(concentrations, alphaX, alphaY, bcLeft, bcRight, bcTop, bcBottom, colMax, i, sx, sy)
|
|
||||||
|
|
||||||
row_t1 = solverFunc(A, b)
|
|
||||||
|
|
||||||
concentrations_t1[i, :] = row_t1
|
|
||||||
end
|
end
|
||||||
|
|
||||||
concentrations_t1 = copy(transpose(concentrations_t1))
|
concentrations_intermediate = copy(transpose(concentrations_intermediate))
|
||||||
concentrations = copy(transpose(concentrations))
|
concentrations_t = fetch(concentrations_t_task)
|
||||||
alphaX = getAlphaX_t(grid)
|
|
||||||
alphaY = getAlphaY_t(grid)
|
|
||||||
|
|
||||||
@inbounds for i = 1:colMax
|
localBs = [zeros(T, rows) for _ in 1:Threads.nthreads()]
|
||||||
|
Threads.@threads for i = 1:cols
|
||||||
|
localB = localBs[Threads.threadid()]
|
||||||
# Swap alphas, boundary conditions and sx/sy for column-wise calculation
|
# Swap alphas, boundary conditions and sx/sy for column-wise calculation
|
||||||
A = createCoeffMatrix(alphaY, bcTop, bcBottom, rowMax, i, sy)
|
A::Tridiagonal{T} = createCoeffMatrix(alphaY_t, alphaY_t_left[i, :], alphaY_t_right[i, :], bcTop, bcBottom, rows, i, sy)
|
||||||
b = createSolutionVector(concentrations_t1, alphaY, alphaX, bcTop, bcBottom, bcLeft, bcRight, rowMax, i, sy, sx)
|
writeSolutionVector!(localB, concentrations_intermediate, alphaY_t, alphaX_t, bcTop, bcBottom, bcLeft, bcRight, i, sy, sx)
|
||||||
|
|
||||||
row_t1 = solverFunc(A, b)
|
concentrations_t[i, :] = (A \ localB)
|
||||||
|
|
||||||
concentrations[i, :] = row_t1
|
|
||||||
end
|
end
|
||||||
|
|
||||||
concentrations = copy(transpose(concentrations))
|
concentrations = copy(transpose(concentrations_t))
|
||||||
|
|
||||||
setConcentrations!(grid, concentrations)
|
setConcentrations!(grid, concentrations)
|
||||||
end
|
end
|
||||||
|
|
||||||
function BTCS_step(grid::Grid{T}, bc::Boundary{T}, timestep::T, numThreads::Int=1) where {T}
|
function runBTCS(grid::Grid{T}, bc::Boundary{T}, timestep::T, iterations::Int, stepCallback::Function) where {T}
|
||||||
if grid.dim == 1
|
if grid.dim == 1
|
||||||
BTCS_1D(grid, bc, timestep, LinearAlgebraAlgorithm)
|
for _ in 1:(iterations)
|
||||||
|
BTCS_1D(grid, bc, timestep)
|
||||||
|
|
||||||
|
stepCallback()
|
||||||
|
end
|
||||||
elseif grid.dim == 2
|
elseif grid.dim == 2
|
||||||
BTCS_2D(grid, bc, timestep, LinearAlgebraAlgorithm, numThreads)
|
alphaX = getAlphaX(grid)
|
||||||
|
alphaY_t = getAlphaY_t(grid)
|
||||||
|
|
||||||
|
alphaX_left_task = Threads.@spawn calcAlphaIntercell(alphaX[:, 1:(grid.cols-2)], alphaX[:, 2:(grid.cols-1)])
|
||||||
|
alphaX_right_task = Threads.@spawn calcAlphaIntercell(alphaX[:, 2:(grid.cols-1)], alphaX[:, 3:grid.cols])
|
||||||
|
alphaY_t_left_task = Threads.@spawn calcAlphaIntercell(alphaY_t[:, 1:(grid.rows-2)], alphaY_t[:, 2:(grid.rows-1)])
|
||||||
|
alphaY_t_right_task = Threads.@spawn calcAlphaIntercell(alphaY_t[:, 2:(grid.rows-1)], alphaY_t[:, 3:grid.rows])
|
||||||
|
|
||||||
|
alphaX_left = fetch(alphaX_left_task)
|
||||||
|
alphaX_right = fetch(alphaX_right_task)
|
||||||
|
alphaY_t_left = fetch(alphaY_t_left_task)
|
||||||
|
alphaY_t_right = fetch(alphaY_t_right_task)
|
||||||
|
|
||||||
|
for _ in 1:(iterations)
|
||||||
|
BTCS_2D(grid, bc, alphaX_left, alphaX_right, alphaY_t_left, alphaY_t_right, timestep)
|
||||||
|
|
||||||
|
stepCallback()
|
||||||
|
end
|
||||||
else
|
else
|
||||||
error("Error: Only 1- and 2-dimensional grids are defined!")
|
error("Error: Only 1- and 2-dimensional grids are defined!")
|
||||||
end
|
end
|
||||||
|
|
||||||
end
|
end
|
||||||
|
|||||||
@ -42,10 +42,6 @@ struct Grid{T}
|
|||||||
end
|
end
|
||||||
end
|
end
|
||||||
|
|
||||||
function setConcentrations!(grid::Grid{T}, new_concentrations::Matrix{T}) where {T}
|
|
||||||
grid.concentrations[] = new_concentrations
|
|
||||||
end
|
|
||||||
|
|
||||||
function getAlphaX(grid::Grid{T})::Matrix{T} where {T}
|
function getAlphaX(grid::Grid{T})::Matrix{T} where {T}
|
||||||
grid.alphaX
|
grid.alphaX
|
||||||
end
|
end
|
||||||
@ -61,3 +57,19 @@ end
|
|||||||
function getAlphaY_t(grid::Grid{T})::Matrix{T} where {T}
|
function getAlphaY_t(grid::Grid{T})::Matrix{T} where {T}
|
||||||
grid.alphaY_t
|
grid.alphaY_t
|
||||||
end
|
end
|
||||||
|
|
||||||
|
function getConcentrations(grid::Grid{T})::Matrix{T} where {T}
|
||||||
|
grid.concentrations[]
|
||||||
|
end
|
||||||
|
|
||||||
|
function setConcentrations!(grid::Grid{T}, new_concentrations::Matrix{T}) where {T}
|
||||||
|
grid.concentrations[] = new_concentrations
|
||||||
|
end
|
||||||
|
|
||||||
|
function getCols(grid::Grid{T})::Int where {T}
|
||||||
|
grid.cols
|
||||||
|
end
|
||||||
|
|
||||||
|
function getRows(grid::Grid{T})::Int where {T}
|
||||||
|
grid.rows
|
||||||
|
end
|
||||||
|
|||||||
@ -90,31 +90,26 @@ function run(simulation::Simulation{T,approach,solver}) where {T,approach,solver
|
|||||||
file = createCSVfile(simulation)
|
file = createCSVfile(simulation)
|
||||||
end
|
end
|
||||||
|
|
||||||
if simulation.approach == BTCS
|
function simulationStepCallback()
|
||||||
if simulation.solver == EIGEN_LU_SOLVER
|
if simulation.consoleOutput >= CONSOLE_OUTPUT_VERBOSE
|
||||||
for _ in 1:(simulation.iterations)
|
printConcentrations(simulation)
|
||||||
if simulation.consoleOutput >= CONSOLE_OUTPUT_VERBOSE
|
|
||||||
printConcentrations(simulation)
|
|
||||||
end
|
|
||||||
|
|
||||||
if simulation.csvOutput >= CSV_OUTPUT_VERBOSE
|
|
||||||
printConcentrationsCSV(simulation, file)
|
|
||||||
end
|
|
||||||
|
|
||||||
BTCS_step(simulation.grid, simulation.bc, simulation.timestep)
|
|
||||||
end
|
|
||||||
else
|
|
||||||
error("Undefined solver!")
|
|
||||||
end
|
end
|
||||||
|
if simulation.csvOutput >= CSV_OUTPUT_VERBOSE
|
||||||
|
printConcentrationsCSV(simulation, file)
|
||||||
|
end
|
||||||
|
end
|
||||||
|
|
||||||
|
if simulation.approach == BTCS
|
||||||
|
runBTCS(simulation.grid, simulation.bc, simulation.timestep, simulation.iterations, simulationStepCallback)
|
||||||
else
|
else
|
||||||
error("Undefined approach!")
|
error("Undefined approach!")
|
||||||
end
|
end
|
||||||
|
|
||||||
if simulation.consoleOutput == CONSOLE_OUTPUT_ON || simulation.consoleOutput == CONSOLE_OUTPUT_VERBOSE
|
if simulation.consoleOutput >= CONSOLE_OUTPUT_ON
|
||||||
printConcentrations(simulation)
|
printConcentrations(simulation)
|
||||||
end
|
end
|
||||||
|
|
||||||
if simulation.csvOutput == CSV_OUTPUT_ON || simulation.csvOutput == CSV_OUTPUT_VERBOSE || simulation.csvOutput == CSV_OUTPUT_XTREME
|
if simulation.csvOutput >= CSV_OUTPUT_ON
|
||||||
printConcentrationsCSV(simulation, file)
|
printConcentrationsCSV(simulation, file)
|
||||||
end
|
end
|
||||||
|
|
||||||
|
|||||||
Loading…
x
Reference in New Issue
Block a user