12#ifndef CD_KMCSolverImplem_H
13#define CD_KMCSolverImplem_H
17#include <unordered_set>
26#include <CD_NamespaceHeader.H>
28template <
typename R,
typename State,
typename T>
31 CH_TIME(
"KMCSolver::KMCSolver");
33 this->setSolverParameters(0, 0, 100, std::numeric_limits<Real>::max(), 0.0, 1.E-6);
36template <
typename R,
typename State,
typename T>
39 CH_TIME(
"KMCSolver::KMCSolver");
44template <
typename R,
typename State,
typename T>
47 CH_TIME(
"KMCSolver::~KMCSolver");
50template <
typename R,
typename State,
typename T>
59 this->setSolverParameters(0, 0, 100, std::numeric_limits<Real>::max(), 0.0, 1.E-6);
62template <
typename R,
typename State,
typename T>
71 CH_TIME(
"KMCSolver::setSolverParameters");
81template <
typename R,
typename State,
typename T>
85 CH_TIME(
"KMCSolver::getNu(State)");
97 x =
static_cast<T
>(0);
116template <
typename R,
typename State,
typename T>
120 CH_TIME(
"KMCSolver::propensities(State)");
122 return this->propensities(
a_state, m_reactions);
125template <
typename R,
typename State,
typename T>
129 CH_TIME(
"KMCSolver::propensities(State, ReactionList)");
142template <
typename R,
typename State,
typename T>
146 CH_TIME(
"KMCSolver::totalPropensity(State)");
148 return this->totalPropensity(
a_state, m_reactions);
151template <
typename R,
typename State,
typename T>
155 CH_TIME(
"KMCSolver::totalPropensity(State, ReactionList)");
168template <
typename R,
typename State,
typename T>
172 CH_TIME(
"KMCSolver::partitionReactions(State)");
174 return this->partitionReactions(
a_state, m_reactions);
177template <
typename R,
typename State,
typename T>
181 CH_TIME(
"KMCSolver::partitionReactions(State, ReactionList)");
202template <
typename R,
typename State,
typename T>
207 return this->getCriticalTimeStep(
a_state, m_reactions);
210template <
typename R,
typename State,
typename T>
215 CH_TIME(
"KMCSolver::getCriticalTimeStep");
219 Real dt = std::numeric_limits<Real>::max();
225 dt = this->getCriticalTimeStep(
A);
231template <
typename R,
typename State,
typename T>
235 CH_TIME(
"KineticMonteCarlo::getCriticalTimeStep(std::vector<Real>)");
238 Real A = std::numeric_limits<Real>::min();
244 return this->getCriticalTimeStep(
A);
247template <
typename R,
typename State,
typename T>
251 CH_TIME(
"KineticMonteCarlo::getCriticalTimeStep(Real)");
258template <
typename R,
typename State,
typename T>
262 CH_TIME(
"KMCSolver::getNonCriticalTimeStep(State");
269template <
typename R,
typename State,
typename T>
273 CH_TIME(
"KMCSolver::getNonCriticalTimeStep(State, ReactionList");
280template <
typename R,
typename State,
typename T>
286 CH_TIME(
"KMCSolver::getNonCriticalTimeStep(State, ReactionList, std::vector<Real>");
292 Real dt = std::numeric_limits<Real>::max();
336 Real dt1 = std::numeric_limits<Real>::max();
337 Real dt2 = std::numeric_limits<Real>::max();
344 if (
mu > std::numeric_limits<Real>::min()) {
347 if (
sigma2 > std::numeric_limits<Real>::min()) {
359template <
typename R,
typename State,
typename T>
366 CH_TIME(
"KMCSolver::computeDt(State, ReactionList, std::vector<Real>, Real");
372 Real dt = std::numeric_limits<Real>::max();
418 if (
mu > std::numeric_limits<Real>::min()) {
428template <
typename R,
typename State,
typename T>
432 CH_TIME(
"KMCSolver::stepSSA(State, Real)");
434 this->stepSSA(
a_state, m_reactions);
437template <
typename R,
typename State,
typename T>
441 CH_TIME(
"KMCSolver::stepSSA(State, ReactionList, Real)");
452template <
typename R,
typename State,
typename T>
458 CH_TIME(
"KMCSolver::stepSSA(State, ReactionList, std::vector<Real>)");
465 constexpr T
one = (T)1;
475 size_t r = std::numeric_limits<size_t>::max();
493template <
typename R,
typename State,
typename T>
497 CH_TIME(
"KMCSolver::advanceSSA(State, Real)");
502template <
typename R,
typename State,
typename T>
506 CH_TIME(
"KMCSolver::advanceSSA(State, ReactionList, Real)");
520 const Real nextDt = this->getCriticalTimeStep(propensities);
532template <
typename R,
typename State,
typename T>
536 CH_TIME(
"KMCSolver::stepExplicitEuler(State, Real)");
538 this->stepExplicitEuler(
a_state, m_reactions,
a_dt);
541template <
typename R,
typename State,
typename T>
547 CH_TIME(
"KMCSolver::stepExplicitEuler(State, ReactionList, std::vector<Real>, Real)");
565template <
typename R,
typename State,
typename T>
569 CH_TIME(
"KMCSolver::stepMidpoint(State&, const Real)");
574template <
typename R,
typename State,
typename T>
578 CH_TIMERS(
"KMCSolver::stepMidpoint(State&, const ReactionList&, const Real)");
579 CH_TIMER(
"KMCSolver::stepTauMidPoint::first_loop",
t1);
580 CH_TIMER(
"KMCSolver::stepTauMidPoint::second_loop",
t2);
610template <
typename R,
typename State,
typename T>
614 CH_TIME(
"KMCSolver::stepPRC(State&, const Real)");
619template <
typename R,
typename State,
typename T>
623 CH_TIMERS(
"KMCSolver::stepPRC(State&, const ReactionList&, const Real)");
624 CH_TIMER(
"KMCSolver::stepPRC::first_loop",
t1);
625 CH_TIMER(
"KMCSolver::stepPRC::second_loop",
t2);
651 const T
nr = (T)Random::getPoisson<long long>(
aj[
i] *
a_dt);
659template <
typename R,
typename State,
typename T>
663 CH_TIME(
"KMCSolver::stepImplicitEuler(State, a_dt)");
665 this->stepImplicitEuler(
a_state, m_reactions,
a_dt);
668template <
typename R,
typename State,
typename T>
674 CH_TIME(
"KMCSolver::stepImplicitEuler(State, ReactionList, a_dt)");
706 for (
int i = 0;
i <
N;
i++) {
712 for (
int j = 0;
j <
M;
j++) {
715 for (
int i = 0;
i <
N;
i++) {
722 auto computeF = [&](
double*
F,
const double*
Xit,
const double*
C) ->
void {
727 for (
int i = 0;
i <
N;
i++) {
736 for (
int i = 0;
i <
N;
i++) {
740 for (
int j = 0;
j <
M;
j++) {
743 for (
int i = 0;
i <
N;
i++) {
755 for (
int i = 0;
i <
N;
i++) {
763 double*
J =
new double[
N *
N];
764 double*
X =
new double[
N];
765 double*
F =
new double[
N];
766 double*
C =
new double[
N];
769 double*
X1 =
new double[
N];
770 double*
X2 =
new double[
N];
771 double*
F1 =
new double[
N];
772 double*
F2 =
new double[
N];
777 int*
IPIV =
new int[
N];
783 for (
int i = 0;
i <
N;
i++) {
791 for (
int k = 0;
k < m_maxIter;
k++) {
795 for (
int j = 0;
j <
N;
j++) {
797 for (
int s = 0;
s <
N;
s++) {
802 X2[
j] += std::max(0.01 *
X[
j], 1.0);
807 for (
int i = 0;
i <
N;
i++) {
821 const std::string err =
"KMCSolver<R, State, T>::stepImplicitEuler -- could not solve A*x = b";
831 for (
int i = 0;
i <
N;
i++) {
848 for (
int i = 0;
i <
N;
i++) {
854 for (
int i = 0;
i <
N;
i++) {
877template <
typename R,
typename State,
typename T>
883 CH_TIME(
"KMCSolver::advanceTau(State, ReactionList, Real, KMCLeapPropagator)");
888template <
typename R,
typename State,
typename T>
895 CH_TIME(
"KMCSolver::advanceTau(State, ReactionList, Real, KMCLeapPropagator)");
915 case KMCLeapPropagator::ExplicitEuler: {
920 case KMCLeapPropagator::Midpoint: {
925 case KMCLeapPropagator::PRC: {
930 case KMCLeapPropagator::ImplicitEuler: {
956template <
typename R,
typename State,
typename T>
962 CH_TIME(
"KMCSolver::advanceHybrid(State, Real, KMCLeapPropagator)");
967template <
typename R,
typename State,
typename T>
974 CH_TIME(
"KMCSolver::advanceHybrid(State, ReactionList, Real, KMCLeapPropagator)");
977 case KMCLeapPropagator::ExplicitEuler: {
979 this->stepExplicitEuler(
s,
r,
dt);
984 case KMCLeapPropagator::Midpoint: {
986 this->stepMidpoint(
s,
r,
dt);
991 case KMCLeapPropagator::PRC: {
993 this->stepPRC(
s,
r,
dt);
998 case KMCLeapPropagator::ImplicitEuler: {
1000 this->stepImplicitEuler(
s,
r,
dt);
1006 MayDay::Error(
"KMCSolver::advanceHybrid - unknown leap propagator requested");
1011template <
typename R,
typename State,
typename T>
1019 CH_TIME(
"KMCSolver::advanceHybrid(State, ReactionList, Real, std::function)");
1021 constexpr T
one = (T)1;
1064 if (
A *
curDt < m_SSAlim) {
1078 for (
size_t i = 0;
i < propensities.size();
i++) {
1079 A += propensities[
i];
1142#include <CD_NamespaceFooter.H>
Class for running Kinetic Monte Carlo functionality.
KMCLeapPropagator
Supported propagators for hybrid tau leaping.
Definition CD_KMCSolver.H:33
Interface to some LaPack routines.
void dgesv_(int *N, int *NRHS, double *A, int *LDA, int *IPIV, double *B, int *LDB, int *INFO)
Interface to LaPack for computing solution to Ax=b.
File containing some useful static methods related to random number generation.
virtual ~KMCSolver() noexcept
Destructor.
Definition CD_KMCSolverImplem.H:45
void advanceHybrid(State &a_state, const Real a_dt, const KMCLeapPropagator &a_leapPropagator=KMCLeapPropagator::ExplicitEuler) const noexcept
Advance using Cao et. al. hybrid algorithm over the input time. This can end up using substepping.
Definition CD_KMCSolverImplem.H:958
void stepPRC(State &a_state, const Real a_dt) const noexcept
Perform one leaping step using the PRC method for ALL reactions over a time step a_dt.
Definition CD_KMCSolverImplem.H:612
void stepMidpoint(State &a_state, const Real a_dt) const noexcept
Perform one leaping step using the midpoint method for ALL reactions over a time step a_dt.
Definition CD_KMCSolverImplem.H:567
Real computeDt(const State &a_state, const ReactionList &a_reactions, const std::vector< Real > &a_propensities, const Real a_epsilon) const noexcept
Compute a time step, using the leap condition on the mean value with eps = 1.
Definition CD_KMCSolverImplem.H:361
void define(const ReactionList &a_reactions) noexcept
Define function. Sets the reactions.
Definition CD_KMCSolverImplem.H:52
void setSolverParameters(const T a_numCrit, const T a_numSSA, const T a_maxIter, const Real a_eps, const Real a_SSAlim, const Real a_exitTol) noexcept
Set solver parameters.
Definition CD_KMCSolverImplem.H:64
KMCSolver() noexcept
Default constructor – must subsequently define the object.
Definition CD_KMCSolverImplem.H:29
std::vector< std::vector< T > > getNu(const State &a_state, const ReactionList &a_reactions) const noexcept
Compute the state vector changes for all reactions.
Definition CD_KMCSolverImplem.H:83
Real getCriticalTimeStep(const State &a_state) const noexcept
Get the time to the next critical reaction.
Definition CD_KMCSolverImplem.H:204
void stepExplicitEuler(State &a_state, const Real a_dt) const noexcept
Perform one plain tau-leaping step using ALL reactions.
Definition CD_KMCSolverImplem.H:534
void advanceSSA(State &a_state, const Real a_dt) const noexcept
Advance with the SSA over the input time. This can end up using substepping.
Definition CD_KMCSolverImplem.H:495
std::pair< ReactionList, ReactionList > partitionReactions(const State &a_state) const noexcept
Partition reactions into critical and non-critical reactions.
Definition CD_KMCSolverImplem.H:170
void stepSSA(State &a_state) const noexcept
Perform a single SSA step.
Definition CD_KMCSolverImplem.H:430
void advanceTau(State &a_state, const Real &a_dt, const KMCLeapPropagator &a_leapPropagator=KMCLeapPropagator::ExplicitEuler) const noexcept
Advance using a specified tau-leaping algorithm.
Definition CD_KMCSolverImplem.H:879
Real getNonCriticalTimeStep(const State &a_state) const noexcept
Get the non-critical time step.
Definition CD_KMCSolverImplem.H:260
Real totalPropensity(const State &a_state) const noexcept
Compute the total propensity for ALL reactions.
Definition CD_KMCSolverImplem.H:144
void stepImplicitEuler(State &a_state, const Real a_dt) const noexcept
Perform one implicit Euler tau-leaping step using ALL reactions.
Definition CD_KMCSolverImplem.H:661
std::vector< Real > propensities(const State &a_state) const noexcept
Compute propensities for ALL reactions.
Definition CD_KMCSolverImplem.H:118
static Real getUniformReal01()
Get a uniform real number on the interval [0,1].
Definition CD_RandomImplem.H:152
Base class for a tracer particle solver. This solver can advance particles in a pre-defined velocity ...
Definition CD_TracerParticleSolver.H:37
TracerParticleSolver()
Default constructor.
Definition CD_TracerParticleSolverImplem.H:25