/** * Generates a new set of values for the S vector that is a) random * and b) obeys the conservation rules. */ void SteadyState::randomizeInitialCondition( const Eref& me ) { #ifdef USE_GSL Id ksolve = Field< Id >::get( stoich_, "ksolve" ); vector< double > nVec = LookupField< unsigned int, vector< double > >::get( ksolve,"nVec", 0 ); int numConsv = total_.size(); recalcTotal( total_, gamma_, &nVec[0] ); // The reorderRows function likes to have an I matrix at the end of // numVarPools_, so we provide space for it, although only its first // column is used for the total vector. gsl_matrix* U = gsl_matrix_calloc ( numConsv, numVarPools_ + numConsv ); for ( int i = 0; i < numConsv; ++i ) { for ( unsigned int j = 0; j < numVarPools_; ++j ) { gsl_matrix_set( U, i, j, gsl_matrix_get( gamma_, i, j ) ); } gsl_matrix_set( U, i, numVarPools_, total_[i] ); } // Do the forward elimination int rank = myGaussianDecomp( U ); assert( rank = numConsv ); vector< double > eliminatedTotal( numConsv, 0.0 ); for ( int i = 0; i < numConsv; ++i ) { eliminatedTotal[i] = gsl_matrix_get( U, i, numVarPools_ ); } // Put Find a vector Y that fits the consv rules. vector< double > y( numVarPools_, 0.0 ); do { fitConservationRules( U, eliminatedTotal, y ); } while ( !checkAboveZero( y ) ); // Sanity check. Try the new vector with the old gamma and tots for ( int i = 0; i < numConsv; ++i ) { double tot = 0.0; for ( unsigned int j = 0; j < numVarPools_; ++j ) { tot += y[j] * gsl_matrix_get( gamma_, i, j ); } assert( fabs( tot - total_[i] ) / tot < EPSILON ); } // Put the new values into S. // cout << endl; for ( unsigned int j = 0; j < numVarPools_; ++j ) { nVec[j] = y[j]; // cout << y[j] << " "; } LookupField< unsigned int, vector< double > >::set( ksolve,"nVec", 0, nVec ); #endif }
void SteadyState::setupSSmatrix() { #ifdef USE_GSL if ( numVarPools_ == 0 || nReacs_ == 0 ) return; int nTot = numVarPools_ + nReacs_; gsl_matrix* N = gsl_matrix_calloc (numVarPools_, nReacs_); if ( LU_ ) { // Clear out old one. gsl_matrix_free( LU_ ); } LU_ = gsl_matrix_calloc (numVarPools_, nTot); vector< int > entry = Field< vector< int > >::get( stoich_, "matrixEntry" ); vector< unsigned int > colIndex = Field< vector< unsigned int > >::get( stoich_, "columnIndex" ); vector< unsigned int > rowStart = Field< vector< unsigned int > >::get( stoich_, "rowStart" ); // cout << endl << endl; for ( unsigned int i = 0; i < numVarPools_; ++i ) { gsl_matrix_set (LU_, i, i + nReacs_, 1 ); unsigned int k = rowStart[i]; // cout << endl << i << ": "; for ( unsigned int j = 0; j < nReacs_; ++j ) { double x = 0; if ( j == colIndex[k] && k < rowStart[i+1] ) { x = entry[k++]; } // cout << " " << x; gsl_matrix_set (N, i, j, x); gsl_matrix_set (LU_, i, j, x ); } } cout << endl << endl; rank_ = myGaussianDecomp( LU_ ); unsigned int nConsv = numVarPools_ - rank_; if ( nConsv == 0 ) { cout << "SteadyState::setupSSmatrix(): Number of conserved species == 0. Aborting\n"; return; } if ( Nr_ ) { // Clear out old one. gsl_matrix_free( Nr_ ); } Nr_ = gsl_matrix_calloc ( rank_, nReacs_ ); // Fill up Nr. for ( unsigned int i = 0; i < rank_; i++) for ( unsigned int j = i; j < nReacs_; j++) gsl_matrix_set (Nr_, i, j, gsl_matrix_get( LU_, i, j ) ); if ( gamma_ ) { // Clear out old one. gsl_matrix_free( gamma_ ); } gamma_ = gsl_matrix_calloc (nConsv, numVarPools_ ); // Fill up gamma for ( unsigned int i = rank_; i < numVarPools_; ++i ) for ( unsigned int j = 0; j < numVarPools_; ++j ) gsl_matrix_set( gamma_, i - rank_, j, gsl_matrix_get( LU_, i, j + nReacs_ ) ); // Fill up boundary condition values total_.resize( nConsv ); total_.assign( nConsv, 0.0 ); /* cout << "S = ("; for ( unsigned int j = 0; j < numVarPools_; ++j ) cout << s_->S()[ j ] << ", "; cout << "), Sinit = ( "; for ( unsigned int j = 0; j < numVarPools_; ++j ) cout << s_->Sinit()[ j ] << ", "; cout << ")\n"; */ Id ksolve = Field< Id >::get( stoich_, "ksolve" ); vector< double > nVec = LookupField< unsigned int, vector< double > >::get( ksolve,"nVec", 0 ); if ( nVec.size() >= numVarPools_ ) { for ( unsigned int i = 0; i < nConsv; ++i ) for ( unsigned int j = 0; j < numVarPools_; ++j ) total_[i] += gsl_matrix_get( gamma_, i, j ) * nVec[ j ]; isSetup_ = 1; } else { cout << "Error: SteadyState::setupSSmatrix(): unable to get" "pool numbers from ksolve.\n"; isSetup_ = 0; } gsl_matrix_free( N ); #endif }