/** * Updates the system according to the probabilistic * number of firings mK[i] of each reaction i */ bool CTauLeapMethod::updateSystem() { // Save the current state in case we need to role back. CVector< C_FLOAT64 > OldState = mContainerState; CMathReaction * pReaction = mReactions.array(); CMathReaction * pReactionEnd = pReaction + mNumReactions; const C_FLOAT64 * pK = mK.array(); const C_FLOAT64 * pKEnd = pK + mNumReactions; for (; pReaction != pReactionEnd; ++pK, ++pReaction) { pReaction->fireMultiple(*pK); } const C_FLOAT64 * pSpecies = mContainerState.array() + mFirstReactionSpeciesIndex; const C_FLOAT64 * pSpeciesEnd = pSpecies + mNumReactionSpecies; for (; pSpecies != pSpeciesEnd; ++pSpecies) { if (*pSpecies < -0.5) { // We need to undo the changes mContainerState = OldState; return false; } } return true; }
/** * Updates the priority queue. * * @param rIndex A size_t giving the index of the fired reaction (-1, if no * stochastic reaction has fired) * @param time A C_FLOAT64 holding the current time */ void CHybridMethod::updatePriorityQueue(size_t rIndex, C_FLOAT64 time) { C_FLOAT64 newTime; // Update all propensities depending on the current reaction (rIndex) if (rIndex != C_INVALID_INDEX) { mpContainer->applyUpdateSequence(mUpdateSequences[rIndex]); } else { mpContainer->updateSimulatedValues(false); CMathReaction * pReaction = mReactions.array(); CMathReaction * pReactionEnd = pReaction + mReactions.size(); for (; pReaction != pReactionEnd; ++pReaction) { const_cast< CMathObject * >(pReaction->getPropensityObject())->calculateValue(); } } std::vector< CHybridStochFlag >::const_iterator it = mReactionFlags.begin(); std::vector< CHybridStochFlag >::const_iterator end = mReactionFlags.end(); for (; it != end; ++it) { if (it->mpPrev == NULL) // reaction is stochastic! { mAmuOld[it->mIndex] = mAmu[it->mIndex]; mAmu[it->mIndex] = * (C_FLOAT64 *) mReactions[it->mIndex].getPropensityObject()->getValuePointer(); if (it->mIndex != rIndex) updateTauMu(it->mIndex, time); } } // draw new random number and update the reaction just fired if (rIndex != C_INVALID_INDEX && mReactionFlags[rIndex].mpPrev == NULL) { // reaction is stochastic newTime = time + generateReactionTime(rIndex); mPQ.updateNode(rIndex, newTime); } return; }
void CStochDirectMethod::start() { CTrajectoryMethod::start(); /* get configuration data */ mMaxSteps = getValue< C_INT32 >("Max Internal Steps"); mpRandomGenerator = &mpContainer->getRandomGenerator(); if (getValue< bool >("Use Random Seed")) { mpRandomGenerator->initialize(getValue< unsigned C_INT32 >("Random Seed")); } //mpCurrentState is initialized. This state is not used internally in the //stochastic solver, but it is used for returning the result after each step. //========Initialize Roots Related Arguments======== mNumRoot = mpContainer->getRoots().size(); mRootsFound.resize(mNumRoot); mRootsA.resize(mNumRoot); mRootsB.resize(mNumRoot); mpRootValueNew = &mRootsA; mpRootValueOld = &mRootsB; mRootsNonZero.resize(mNumRoot); mRootsNonZero = 0.0; CMathObject * pRootObject = mpContainer->getMathObject(mpContainer->getRoots().array()); CMathObject * pRootObjectEnd = pRootObject + mNumRoot; CObjectInterface::ObjectSet Requested; for (; pRootObject != pRootObjectEnd; ++pRootObject) { Requested.insert(pRootObject); } CObjectInterface::ObjectSet Changed; // Determine whether we have time dependent roots; CObjectInterface * pTimeObject = mpContainer->getMathObject(mpContainerStateTime); Changed.insert(pTimeObject); mpContainer->getTransientDependencies().getUpdateSequence(mUpdateTimeDependentRoots, CMath::Default, Changed, Requested); mHaveTimeDependentRoots = (mUpdateTimeDependentRoots.size() > 0); // Build the reaction dependencies mReactions.initialize(mpContainer->getReactions()); mNumReactions = mReactions.size(); mAmu.initialize(mpContainer->getPropensities()); mPropensityObjects.initialize(mAmu.size(), mpContainer->getMathObject(mAmu.array())); mUpdateSequences.resize(mNumReactions); CMathReaction * pReaction = mReactions.array(); CMathReaction * pReactionEnd = pReaction + mNumReactions; CObjectInterface::UpdateSequence * pUpdateSequence = mUpdateSequences.array(); CMathObject * pPropensityObject = mPropensityObjects.array(); CMathObject * pPropensityObjectEnd = pPropensityObject + mPropensityObjects.size(); for (; pPropensityObject != pPropensityObjectEnd; ++pPropensityObject) { Requested.insert(pPropensityObject); } pPropensityObject = mPropensityObjects.array(); for (; pReaction != pReactionEnd; ++pReaction, ++pUpdateSequence, ++pPropensityObject) { Changed = pReaction->getChangedObjects(); // The time is always updated Changed.insert(pTimeObject); pUpdateSequence->clear(); mpContainer->getTransientDependencies().getUpdateSequence(*pUpdateSequence, CMath::Default, Changed, Requested); } mMaxStepsReached = false; mTargetTime = *mpContainerStateTime; mNextReactionTime = *mpContainerStateTime; mNextReactionIndex = C_INVALID_INDEX; stateChange(CMath::State); return; }
// virtual void CTrajectoryMethodDsaLsodar::start() { CLsodaMethod::start(); mReactions.initialize(mpContainer->getReactions()); mNumReactions = mReactions.size(); mAmu.initialize(mpContainer->getPropensities()); mPropensityObjects.initialize(mNumReactions, mpContainer->getMathObject(mAmu.array())); mUpdateSequences.resize(mNumReactions); mFirstReactionSpeciesIndex = mpContainer->getCountFixedEventTargets() + 1 /* Time */ + mpContainer->getCountODEs(); // Create a local copy of the state where the particle number species determined // by reactions are rounded to integers. C_FLOAT64 * pValue = mContainerState.array() + mFirstReactionSpeciesIndex; C_FLOAT64 * pValueEnd = pValue + mpContainer->getCountIndependentSpecies() + mpContainer->getCountDependentSpecies(); for (; pValue != pValueEnd; ++pValue) { *pValue = floor(*pValue + 0.5); } // The container state is now up to date we just need to calculate all values needed for simulation. mpContainer->updateSimulatedValues(false); CMathObject * pTimeObject = mpContainer->getMathObject(mpContainer->getModel().getValueReference()); // Build the reaction dependencies mReactions.initialize(mpContainer->getReactions()); mNumReactions = mReactions.size(); mAmu.initialize(mpContainer->getPropensities()); mPropensityObjects.initialize(mAmu.size(), mpContainer->getMathObject(mAmu.array())); mUpdateSequences.resize(mNumReactions); C_FLOAT64 * pAmu = mAmu.array(); mA0 = 0.0; CMathReaction * pReaction = mReactions.array(); CMathReaction * pReactionEnd = pReaction + mNumReactions; CCore::CUpdateSequence * pUpdateSequence; CMathObject * pPropensityObject = mPropensityObjects.array(); CMathObject * pPropensityObjectEnd = pPropensityObject + mPropensityObjects.size(); CObjectInterface::ObjectSet Requested; for (; pPropensityObject != pPropensityObjectEnd; ++pPropensityObject) { Requested.insert(pPropensityObject); } pPropensityObject = mPropensityObjects.array(); for (; pReaction != pReactionEnd; ++pReaction, ++pUpdateSequence, ++pPropensityObject, ++pAmu) { // Update the propensity pPropensityObject->calculateValue(); mA0 += *pAmu; CObjectInterface::ObjectSet Changed; // The time is always updated Changed.insert(pTimeObject); const CMathReaction::SpeciesBalance * itBalance = pReaction->getNumberBalance().array(); const CMathReaction::SpeciesBalance * endBalance = itBalance + pReaction->getNumberBalance().size(); for (; itBalance != endBalance; ++itBalance) { Changed.insert(mpContainer->getMathObject(itBalance->first)); } pUpdateSequence->clear(); mpContainer->getTransientDependencies().getUpdateSequence(*pUpdateSequence, CCore::SimulationContext::Default, Changed, Requested); } mPartition.intialize(mpContainer, *mpLowerLimit, *mpUpperLimit); return; }
/** * Simulates the system over the next interval of time. The timestep * is given as argument. * * @param ds A C_FLOAT64 specifying the timestep */ C_FLOAT64 CTauLeapMethod::doSingleStep(C_FLOAT64 ds) { C_FLOAT64 Lambda, Tmp, Tau, Tau1, Tau2; updatePropensities(); mAvgDX = 0.0; mSigDX = 0.0; CMathReaction * pReaction = mReactions.array(); const C_FLOAT64 * pAmu = mAmu.array(); const C_FLOAT64 * pAmuEnd = pAmu + mNumReactions; const C_FLOAT64 * pFirstSpecies = mContainerState.array() + mFirstReactionSpeciesIndex; for (; pAmu != pAmuEnd; ++pAmu, ++pReaction) { const CMathReaction::Balance & Balance = pReaction->getNumberBalance(); const CMathReaction::SpeciesBalance * it = Balance.array(); const CMathReaction::SpeciesBalance * end = it + Balance.size(); for (; it != end; ++it) { mAvgDX[it->first - pFirstSpecies] += it->second **pAmu; mSigDX[it->first - pFirstSpecies] += it->second * it->second **pAmu; } } Tau1 = Tau2 = std::numeric_limits< C_FLOAT64 >::infinity(); const C_FLOAT64 * pSpecies = mContainerState.array() + mFirstReactionSpeciesIndex; const C_FLOAT64 * pSpeciesEnd = pSpecies + mNumReactionSpecies; C_FLOAT64 * pAvgDX = mAvgDX.array(); C_FLOAT64 * pSigDX = mSigDX.array(); for (; pSpecies != pSpeciesEnd; ++pSpecies, ++pAvgDX, ++pSigDX) { if ((Tmp = mEpsilon * fabs(*pSpecies)) < 1.0) Tmp = 1.0; *pAvgDX = Tmp / fabs(*pAvgDX); *pSigDX = (Tmp * Tmp) / fabs(*pSigDX); if (Tau1 > *pAvgDX) Tau1 = *pAvgDX; if (Tau2 > *pSigDX) Tau2 = *pSigDX; } Tau = std::min(Tau1, Tau2); if (ds < Tau) Tau = ds; pAmu = mAmu.array(); C_FLOAT64 * pK = mK.array(); C_FLOAT64 * pKEnd = pK + mNumReactions; for (; pAmu != pAmuEnd; ++pAmu, ++pK) { Lambda = *pAmu * Tau; if (Lambda < 0.0) CCopasiMessage(CCopasiMessage::EXCEPTION, MCTrajectoryMethod + 10); else if (Lambda > 2.0e9) CCopasiMessage(CCopasiMessage::EXCEPTION, MCTrajectoryMethod + 26); *pK = mpRandomGenerator->getRandomPoisson(Lambda); } while (!updateSystem()) { Tau *= 0.5; pK = mK.array(); for (; pK != pKEnd; ++pK) { *pK *= 0.5; if (*pK < floor(*pK + 0.75)) { *pK += mpRandomGenerator->getRandomCC() < 0.5 ? - 0.5 : 0.5; } } } return Tau; }