Skip to content

Commit

Permalink
195 backward euler python (#234)
Browse files Browse the repository at this point in the history
* Updated python wrapper to include Backward Euler options.

* Added TestAnalyticalBackwardEuler.

* Added TestAnalyticalStandardBackwardEuler.

* Added multiple grid cell backward Euler tests.

* Added places argument with default of 5.

* Fixed wrapper.

---------

Co-authored-by: David Fillmore <fillmore.david.winslow@gmail.com>
  • Loading branch information
dwfncar and davidfillmore authored Nov 5, 2024
1 parent 5b2e2f5 commit 5817b5a
Show file tree
Hide file tree
Showing 2 changed files with 72 additions and 18 deletions.
68 changes: 52 additions & 16 deletions python/test/test_analytical.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,7 @@
import random


def TestSingleGridCell(self, solver):
def TestSingleGridCell(self, solver, places=5):
num_grid_cells = 1
time_step = 200.0
temperature = 272.5
Expand Down Expand Up @@ -78,17 +78,17 @@ def TestSingleGridCell(self, solver):
(1.0 + (k1 * np.exp(-k2 * curr_time) - k2 * np.exp(-k1 * curr_time)) / (k2 - k1))

self.assertAlmostEqual(
ordered_concentrations[species_ordering["A"]], A_conc, places=5)
ordered_concentrations[species_ordering["A"]], A_conc, places=places)
self.assertAlmostEqual(
ordered_concentrations[species_ordering["B"]], B_conc, places=5)
ordered_concentrations[species_ordering["B"]], B_conc, places=places)
self.assertAlmostEqual(
ordered_concentrations[species_ordering["C"]], C_conc, places=5)
ordered_concentrations[species_ordering["C"]], C_conc, places=places)
self.assertAlmostEqual(
ordered_concentrations[species_ordering["D"]], D_conc, places=5)
ordered_concentrations[species_ordering["D"]], D_conc, places=places)
self.assertAlmostEqual(
ordered_concentrations[species_ordering["E"]], E_conc, places=5)
ordered_concentrations[species_ordering["E"]], E_conc, places=places)
self.assertAlmostEqual(
ordered_concentrations[species_ordering["F"]], F_conc, places=5)
ordered_concentrations[species_ordering["F"]], F_conc, places=places)

curr_time += time_step

Expand All @@ -102,7 +102,7 @@ def test_simulation(self):
TestSingleGridCell(self, solver)


class TestAnalyicalStandardRosenbrock(unittest.TestCase):
class TestAnalyticalStandardRosenbrock(unittest.TestCase):
def test_simulation(self):
solver = musica.create_solver(
"configs/analytical",
Expand All @@ -111,7 +111,25 @@ def test_simulation(self):
TestSingleGridCell(self, solver)


def TestMultipleGridCell(self, solver, num_grid_cells):
class TestAnalyticalBackwardEuler(unittest.TestCase):
def test_simulation(self):
solver = musica.create_solver(
"configs/analytical",
musica.micmsolver.backward_euler,
1)
TestSingleGridCell(self, solver, places=2)


class TestAnalyticalStandardBackwardEuler(unittest.TestCase):
def test_simulation(self):
solver = musica.create_solver(
"configs/analytical",
musica.micmsolver.backward_euler_standard_order,
1)
TestSingleGridCell(self, solver, places=2)


def TestMultipleGridCell(self, solver, num_grid_cells, places=5):
time_step = 200.0
temperatures = []
pressures = []
Expand Down Expand Up @@ -221,17 +239,17 @@ def TestMultipleGridCell(self, solver, num_grid_cells):
k1[i] * np.exp(-k2[i] * curr_time) - k2[i] * np.exp(-k1[i] * curr_time)) / (k2[i] - k1[i]))

self.assertAlmostEqual(
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["A"]], A_conc, places=5)
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["A"]], A_conc, places=places)
self.assertAlmostEqual(
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["B"]], B_conc, places=5)
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["B"]], B_conc, places=places)
self.assertAlmostEqual(
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["C"]], C_conc, places=5)
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["C"]], C_conc, places=places)
self.assertAlmostEqual(
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["D"]], D_conc, places=5)
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["D"]], D_conc, places=places)
self.assertAlmostEqual(
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["E"]], E_conc, places=5)
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["E"]], E_conc, places=places)
self.assertAlmostEqual(
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["F"]], F_conc, places=5)
ordered_concentrations[i * len(concentrations.keys()) + species_ordering["F"]], F_conc, places=places)

curr_time += time_step

Expand All @@ -245,7 +263,7 @@ def test_simulation(self):
TestMultipleGridCell(self, solver, 3)


class TestAnalyicalStandardRosenbrockMultipleGridCells(unittest.TestCase):
class TestAnalyticalStandardRosenbrockMultipleGridCells(unittest.TestCase):
def test_simulation(self):
solver = musica.create_solver(
"configs/analytical",
Expand All @@ -254,5 +272,23 @@ def test_simulation(self):
TestMultipleGridCell(self, solver, 3)


class TestAnalyticalBackwardEulerMultipleGridCells(unittest.TestCase):
def test_simulation(self):
solver = musica.create_solver(
"configs/analytical",
musica.micmsolver.backward_euler,
3)
TestMultipleGridCell(self, solver, 3, places=2)


class TestAnalyticalStandardBackwardEulerMultipleGridCells(unittest.TestCase):
def test_simulation(self):
solver = musica.create_solver(
"configs/analytical",
musica.micmsolver.backward_euler_standard_order,
3)
TestMultipleGridCell(self, solver, 3, places=2)


if __name__ == '__main__':
unittest.main()
22 changes: 20 additions & 2 deletions python/wrapper.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -14,7 +14,9 @@ PYBIND11_MODULE(musica, m)

py::enum_<musica::MICMSolver>(m, "micmsolver")
.value("rosenbrock", musica::MICMSolver::Rosenbrock)
.value("rosenbrock_standard_order", musica::MICMSolver::RosenbrockStandardOrder);
.value("rosenbrock_standard_order", musica::MICMSolver::RosenbrockStandardOrder)
.value("backward_euler", musica::MICMSolver::BackwardEuler)
.value("backward_euler_standard_order", musica::MICMSolver::BackwardEulerStandardOrder);

m.def(
"create_solver",
Expand Down Expand Up @@ -162,6 +164,14 @@ PYBIND11_MODULE(musica, m)
{
map = micm->rosenbrock_standard_.second.variable_map_;
}
else if (micm->solver_type_ == musica::MICMSolver::BackwardEuler)
{
map = micm->backward_euler_.second.variable_map_;
}
else if (micm->solver_type_ == musica::MICMSolver::BackwardEulerStandardOrder)
{
map = micm->backward_euler_standard_.second.variable_map_;
}

return map;
},
Expand All @@ -182,8 +192,16 @@ PYBIND11_MODULE(musica, m)
{
map = micm->rosenbrock_standard_.second.custom_rate_parameter_map_;
}
else if (micm->solver_type_ == musica::MICMSolver::BackwardEuler)
{
map = micm->backward_euler_.second.custom_rate_parameter_map_;
}
else if (micm->solver_type_ == musica::MICMSolver::BackwardEulerStandardOrder)
{
map = micm->backward_euler_standard_.second.custom_rate_parameter_map_;
}

return map;
},
"Return map of reaction rates");
}
}

0 comments on commit 5817b5a

Please sign in to comment.