-
Notifications
You must be signed in to change notification settings - Fork 3
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
196 backward euler fortran #210
Merged
Merged
Changes from 10 commits
Commits
Show all changes
12 commits
Select commit
Hold shift + click to select a range
9fb15b6
Added Backward Euler definitions to musics/mich.hpp.
dwfncar 33e55a4
Added Backward Euler definitions micm.cpp.
dwfncar dea6f4d
Added unit test BackwardEulerStandardOrder.
dwfncar 5f08aec
Added unit test BackwardEulerVectorOrder.
dwfncar c6fcf7a
Save.
dwfncar efa5d8e
Added Backward Euler fortran API.
dwfncar df2a1e3
Merged main.
davidfillmore 414a14e
Added accuracy parameter for ASSERT_NEAR.
dwfncar fbd9581
Reduced time_step for Backward Euler.
dwfncar 16bfed5
Added time_step and test_accuracy arguments to test_mulitple_grid_cells.
dwfncar e8b8c5e
Removed some prints from test_micm_api.F90.
dwfncar d891a9b
Merged main to 196-backward-euler-fortran.
dwfncar 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 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 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 |
---|---|---|
|
@@ -6,7 +6,7 @@ program test_micm_api | |
use, intrinsic :: iso_c_binding | ||
use, intrinsic :: ieee_arithmetic | ||
use musica_micm, only: micm_t, solver_stats_t, get_micm_version | ||
use musica_micm, only: Rosenbrock, RosenbrockStandardOrder | ||
use musica_micm, only: Rosenbrock, RosenbrockStandardOrder, BackwardEuler, BackwardEulerStandardOrder | ||
use musica_util, only: assert, error_t, mapping_t, string_t, find_mapping_index | ||
|
||
#include "micm/util/error.hpp" | ||
|
@@ -29,6 +29,8 @@ program test_micm_api | |
call test_api() | ||
call test_multiple_grid_cell_vector_Rosenbrock() | ||
call test_multiple_grid_cell_standard_Rosenbrock() | ||
call test_multiple_grid_cell_vector_BackwardEuler() | ||
call test_multiple_grid_cell_standard_BackwardEuler() | ||
|
||
contains | ||
|
||
|
@@ -170,14 +172,15 @@ subroutine test_api() | |
|
||
end subroutine test_api | ||
|
||
subroutine test_multiple_grid_cells(micm, NUM_GRID_CELLS) | ||
subroutine test_multiple_grid_cells(micm, NUM_GRID_CELLS, time_step, test_accuracy) | ||
|
||
type(micm_t), pointer, intent(inout) :: micm | ||
integer, intent(in) :: NUM_GRID_CELLS | ||
real(c_double), intent(in) :: time_step | ||
real, intent(in) :: test_accuracy | ||
|
||
integer, parameter :: NUM_SPECIES = 6 | ||
integer, parameter :: NUM_USER_DEFINED_REACTION_RATES = 2 | ||
real(c_double) :: time_step | ||
real(c_double), target :: temperature(NUM_GRID_CELLS) | ||
real(c_double), target :: pressure(NUM_GRID_CELLS) | ||
real(c_double), target :: air_density(NUM_GRID_CELLS) | ||
|
@@ -199,8 +202,6 @@ subroutine test_multiple_grid_cells(micm, NUM_GRID_CELLS) | |
real :: temp | ||
type(ArrheniusReaction) :: r1, r2 | ||
|
||
time_step = 200 | ||
|
||
A_index = find_mapping_index( micm%species_ordering, "A", found ) | ||
ASSERT( found ) | ||
B_index = find_mapping_index( micm%species_ordering, "B", found ) | ||
|
@@ -271,12 +272,19 @@ subroutine test_multiple_grid_cells(micm, NUM_GRID_CELLS) | |
D = initial_D * exp( -k1 * time_step ) | ||
E = initial_D * (k1 / (k2 - k1)) * (exp(-k1 * time_step) - exp(-k2 * time_step)) | ||
F = initial_F + initial_D * (1.0 + (k1 * exp(-k2 * time_step) - k2 * exp(-k1 * time_step)) / (k2 - k1)) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+A_index), A, 5.0e-3) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+B_index), B, 5.0e-3) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+C_index), C, 5.0e-3) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+D_index), D, 5.0e-3) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+E_index), E, 5.0e-3) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+F_index), F, 5.0e-3) | ||
print *, 'A', concentrations((i_cell-1)*NUM_SPECIES+A_index), A | ||
print *, 'B', concentrations((i_cell-1)*NUM_SPECIES+B_index), B | ||
print *, 'C', concentrations((i_cell-1)*NUM_SPECIES+C_index), C | ||
print *, 'D', concentrations((i_cell-1)*NUM_SPECIES+D_index), D | ||
print *, 'E', concentrations((i_cell-1)*NUM_SPECIES+E_index), E | ||
print *, 'F', concentrations((i_cell-1)*NUM_SPECIES+F_index), F | ||
print * | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+A_index), A, test_accuracy) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+B_index), B, test_accuracy) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+C_index), C, test_accuracy) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+D_index), D, test_accuracy) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+E_index), E, test_accuracy) | ||
ASSERT_NEAR(concentrations((i_cell-1)*NUM_SPECIES+F_index), F, test_accuracy) | ||
end do | ||
|
||
end subroutine test_multiple_grid_cells | ||
|
@@ -287,11 +295,15 @@ subroutine test_multiple_grid_cell_vector_Rosenbrock() | |
integer :: num_grid_cells | ||
type(error_t) :: error | ||
|
||
real(c_double), parameter :: time_step = 200 | ||
real, parameter :: test_accuracy = 5.0e-3 | ||
|
||
num_grid_cells = 3 | ||
micm => micm_t( "configs/analytical", Rosenbrock, num_grid_cells, error ) | ||
ASSERT( error%is_success() ) | ||
|
||
call test_multiple_grid_cells( micm, num_grid_cells ) | ||
print *, 'test_multiple_grid_cells vector Rosenbrock' | ||
call test_multiple_grid_cells( micm, num_grid_cells, time_step, test_accuracy ) | ||
|
||
deallocate( micm ) | ||
|
||
|
@@ -303,14 +315,58 @@ subroutine test_multiple_grid_cell_standard_Rosenbrock() | |
integer :: num_grid_cells | ||
type(error_t) :: error | ||
|
||
real(c_double), parameter :: time_step = 200 | ||
real, parameter :: test_accuracy = 5.0e-3 | ||
|
||
num_grid_cells = 3 | ||
micm => micm_t( "configs/analytical", RosenbrockStandardOrder, num_grid_cells, error ) | ||
ASSERT( error%is_success() ) | ||
|
||
call test_multiple_grid_cells( micm, num_grid_cells ) | ||
print *, 'test_multiple_grid_cells standard Rosenbrock' | ||
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. please remove diagnostic output |
||
call test_multiple_grid_cells( micm, num_grid_cells, time_step, test_accuracy ) | ||
|
||
deallocate( micm ) | ||
|
||
end subroutine test_multiple_grid_cell_standard_Rosenbrock | ||
|
||
subroutine test_multiple_grid_cell_vector_BackwardEuler() | ||
|
||
type(micm_t), pointer :: micm | ||
integer :: num_grid_cells | ||
type(error_t) :: error | ||
|
||
real(c_double), parameter :: time_step = 10 | ||
real, parameter :: test_accuracy = 0.1 | ||
|
||
num_grid_cells = 3 | ||
micm => micm_t( "configs/analytical", BackwardEuler, num_grid_cells, error ) | ||
ASSERT( error%is_success() ) | ||
|
||
print *, 'test_multiple_grid_cells vector Backward Euler' | ||
call test_multiple_grid_cells( micm, num_grid_cells, time_step, test_accuracy ) | ||
|
||
deallocate( micm ) | ||
|
||
end subroutine test_multiple_grid_cell_vector_BackwardEuler | ||
|
||
subroutine test_multiple_grid_cell_standard_BackwardEuler() | ||
|
||
type(micm_t), pointer :: micm | ||
integer :: num_grid_cells | ||
type(error_t) :: error | ||
|
||
real(c_double), parameter :: time_step = 10 | ||
real, parameter :: test_accuracy = 0.1 | ||
|
||
num_grid_cells = 3 | ||
micm => micm_t( "configs/analytical", BackwardEulerStandardOrder, num_grid_cells, error ) | ||
ASSERT( error%is_success() ) | ||
|
||
print *, 'test_multiple_grid_cells standard Backward Euler' | ||
call test_multiple_grid_cells( micm, num_grid_cells, time_step, test_accuracy ) | ||
|
||
deallocate( micm ) | ||
|
||
end subroutine test_multiple_grid_cell_standard_BackwardEuler | ||
|
||
end program |
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.
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.
please remove diagnostic output