void ZombieFunction::innerSetExpr( const Eref& e, string v ) { Function::innerSetExpr( e, v ); if ( _stoich ) { Stoich* s = reinterpret_cast< Stoich* >( _stoich ); s->setFunctionExpr( e, v ); } else { // I had this warning here but it is triggered needlessly when we // do an assignment of the Zombie function. So removed. // cout << "Warning: ZombieFunction::setExpr: Stoich not set.\n"; } }
void SteadyState::setStoich( Id value ) { if ( !value.element()->cinfo()->isA( "Stoich" ) ) { cout << "Error: SteadyState::setStoich: Must be of Stoich class\n"; return; } stoich_ = value; Stoich* stoichPtr = reinterpret_cast< Stoich* >( value.eref().data()); numVarPools_ = Field< unsigned int >::get( stoich_, "numVarPools" ); nReacs_ = Field< unsigned int >::get( stoich_, "numRates" ); setupSSmatrix(); double vol = LookupField< unsigned int, double >::get( stoichPtr->getCompartment(), "oneVoxelVolume", 0 ); pool_.setVolume( vol ); pool_.setStoich( stoichPtr, 0 ); pool_.updateAllRateTerms( stoichPtr->getRateTerms(), stoichPtr->getNumCoreRates() ); isInitialized_ = 1; }
int ss_func( const gsl_vector* x, void* params, gsl_vector* f ) { struct reac_info* ri = (struct reac_info *)params; Stoich* s = reinterpret_cast< Stoich* >( ri->stoich.eref().data() ); int num_consv = ri->num_mols - ri->rank; for ( unsigned int i = 0; i < ri->num_mols; ++i ) { double temp = op( gsl_vector_get( x, i ) ); if ( isnan( temp ) || isinf( temp ) ) { return GSL_ERANGE; } else { ri->nVec[i] = temp; } } vector< double > vels; s->updateReacVelocities( &ri->nVec[0], vels, 0 ); // use compt zero assert( vels.size() == static_cast< unsigned int >( ri->num_reacs ) ); // y = Nr . v // Note that Nr is row-echelon: diagonal and above. for ( int i = 0; i < ri->rank; ++i ) { double temp = 0; for ( int j = i; j < ri->num_reacs; ++j ) temp += gsl_matrix_get( ri->Nr, i, j ) * vels[j]; gsl_vector_set( f, i, temp ); } // dT = gamma.S - T for ( int i = 0; i < num_consv; ++i ) { double dT = - ri->T[i]; for ( unsigned int j = 0; j < ri->num_mols; ++j ) dT += gsl_matrix_get( ri->gamma, i, j) * op( gsl_vector_get( x, j ) ); gsl_vector_set( f, i + ri->rank, dT ); } return GSL_SUCCESS; }
void SteadyState::classifyState( const double* T ) { #ifdef USE_GSL // unsigned int nConsv = numVarPools_ - rank_; gsl_matrix* J = gsl_matrix_calloc ( numVarPools_, numVarPools_ ); // double* yprime = new double[ numVarPools_ ]; // vector< double > yprime( numVarPools_, 0.0 ); // Generate an approximation to the Jacobean by generating small // increments to each of the molecules in the steady state, one // at a time, and putting the resultant rate vector into a column // of the J matrix. // This needs a bit of heuristic to decide what is a 'small' increment. // Use the CoInits for this. Stoichiometry shouldn't matter too much. // I used the totals from consv rules earlier, but that can have // negative values. double tot = 0.0; Stoich* s = reinterpret_cast< Stoich* >( stoich_.eref().data() ); vector< double > nVec = LookupField< unsigned int, vector< double > >::get( s->getKsolve(), "nVec", 0 ); for ( unsigned int i = 0; i < numVarPools_; ++i ) { tot += nVec[i]; } tot *= DELTA; vector< double > yprime( nVec.size(), 0.0 ); // Fill up Jacobian for ( unsigned int i = 0; i < numVarPools_; ++i ) { double orig = nVec[i]; if ( isnan( orig ) ) { cout << "Warning: SteadyState::classifyState: orig=nan\n"; solutionStatus_ = 2; // Steady state OK, eig failed gsl_matrix_free ( J ); return; } if ( isnan( tot ) ) { cout << "Warning: SteadyState::classifyState: tot=nan\n"; solutionStatus_ = 2; // Steady state OK, eig failed gsl_matrix_free ( J ); return; } nVec[i] = orig + tot; /// Here we assume we always use voxel zero. s->updateRates( &nVec[0], &yprime[0], 0 ); nVec[i] = orig; // Assign the rates for each mol. for ( unsigned int j = 0; j < numVarPools_; ++j ) { gsl_matrix_set( J, i, j, yprime[j] ); } } // Jacobian is now ready. Find eigenvalues. gsl_vector_complex* vec = gsl_vector_complex_alloc( numVarPools_ ); gsl_eigen_nonsymm_workspace* workspace = gsl_eigen_nonsymm_alloc( numVarPools_ ); int status = gsl_eigen_nonsymm( J, vec, workspace ); eigenvalues_.clear(); eigenvalues_.resize( numVarPools_, 0.0 ); if ( status != GSL_SUCCESS ) { cout << "Warning: SteadyState::classifyState failed to find eigenvalues. Status = " << status << endl; solutionStatus_ = 2; // Steady state OK, eig classification failed } else { // Eigenvalues are ready. Classify state. nNegEigenvalues_ = 0; nPosEigenvalues_ = 0; for ( unsigned int i = 0; i < numVarPools_; ++i ) { gsl_complex z = gsl_vector_complex_get( vec, i ); double r = GSL_REAL( z ); nNegEigenvalues_ += ( r < -EPSILON ); nPosEigenvalues_ += ( r > EPSILON ); eigenvalues_[i] = r; // We have a problem here because numVarPools_ usually > rank // This means we have several zero eigenvalues. } if ( nNegEigenvalues_ == rank_ ) stateType_ = 0; // Stable else if ( nPosEigenvalues_ == rank_ ) // Never see it. stateType_ = 1; // Unstable else if (nPosEigenvalues_ == 1) stateType_ = 2; // Saddle else if ( nPosEigenvalues_ >= 2 ) stateType_ = 3; // putative oscillatory else if ( nNegEigenvalues_ == ( rank_ - 1) && nPosEigenvalues_ == 0 ) stateType_ = 4; // one zero or unclassified eigenvalue. Messy. else stateType_ = 5; // Other } gsl_vector_complex_free( vec ); gsl_matrix_free ( J ); gsl_eigen_nonsymm_free( workspace ); #endif }
void SteadyState::classifyState( const double* T ) { /* column_major trait is needed for fortran */ ublas::matrix<double, ublas::column_major> J(numVarPools_, numVarPools_); double tot = 0.0; Stoich* s = reinterpret_cast< Stoich* >( stoich_.eref().data() ); vector< double > nVec = LookupField< unsigned int, vector< double > >::get( s->getKsolve(), "nVec", 0 ); for ( unsigned int i = 0; i < numVarPools_; ++i ) tot += nVec[i]; tot *= DELTA; vector< double > yprime( nVec.size(), 0.0 ); // Fill up Jacobian for ( unsigned int i = 0; i < numVarPools_; ++i ) { double orig = nVec[i]; if ( std::isnan( orig ) or std::isinf( orig ) ) { cout << "Warning: SteadyState::classifyState: orig=nan\n"; solutionStatus_ = 2; // Steady state OK, eig failed J.clear(); return; } if ( std::isnan( tot ) or std::isinf( tot )) { cout << "Warning: SteadyState::classifyState: tot=nan\n"; solutionStatus_ = 2; // Steady state OK, eig failed J.clear(); return; } nVec[i] = orig + tot; pool_.updateRates( &nVec[0], &yprime[0] ); nVec[i] = orig; // Assign the rates for each mol. for ( unsigned int j = 0; j < numVarPools_; ++j ) { if( std::isnan( yprime[j] ) or std::isinf( yprime[j] ) ) { cout << "Warning: Overflow/underflow. Can't continue " << endl; solutionStatus_ = 2; return; } J(i, j) = yprime[j]; } } // Jacobian is now ready. Find eigenvalues. ublas::vector< std::complex< double > > eigenVec ( J.size1() ); ublas::matrix< std::complex<double>, ublas::column_major >* vl, *vr; vl = NULL; vr = NULL; /*----------------------------------------------------------------------------- * INFO: Calling lapack routine geev to compute eigen vector of matrix J. * * Argument 3 and 4 are left- and right-eigenvectors. Since we do not need * them, they are set to NULL. Argument 2 holds eigen-vector and result is * copied to it (output ). *-----------------------------------------------------------------------------*/ int status = lapack::geev( J, eigenVec, vl, vr, lapack::optimal_workspace() ); eigenvalues_.clear(); eigenvalues_.resize( numVarPools_, 0.0 ); if ( status != 0 ) { cout << "Warning: SteadyState::classifyState failed to find eigenvalues. Status = " << status << endl; solutionStatus_ = 2; // Steady state OK, eig classification failed } else // Eigenvalues are ready. Classify state. { nNegEigenvalues_ = 0; nPosEigenvalues_ = 0; for ( unsigned int i = 0; i < numVarPools_; ++i ) { std::complex<value_type> z = eigenVec[ i ]; double r = z.real(); nNegEigenvalues_ += ( r < -EPSILON ); nPosEigenvalues_ += ( r > EPSILON ); eigenvalues_[i] = r; // We have a problem here because numVarPools_ usually > rank // This means we have several zero eigenvalues. } if ( nNegEigenvalues_ == rank_ ) stateType_ = 0; // Stable else if ( nPosEigenvalues_ == rank_ ) // Never see it. stateType_ = 1; // Unstable else if (nPosEigenvalues_ == 1) stateType_ = 2; // Saddle else if ( nPosEigenvalues_ >= 2 ) stateType_ = 3; // putative oscillatory else if ( nNegEigenvalues_ == ( rank_ - 1) && nPosEigenvalues_ == 0 ) stateType_ = 4; // one zero or unclassified eigenvalue. Messy. else stateType_ = 5; // Other } }