Browse Source

more commentations, more flexibility for classifier interfaces

Alexander Freytag 11 years ago
parent
commit
bc6844c3a3

+ 50 - 0
FMKGPHyperparameterOptimization.cpp

@@ -1425,6 +1425,56 @@ int FMKGPHyperparameterOptimization::classify ( const NICE::SparseVector & xstar
   }
 }
 
+int FMKGPHyperparameterOptimization::classify ( const NICE::Vector & xstar, NICE::SparseVector & scores ) const
+{
+  // loop through all classes
+  if ( precomputedA.size() == 0 )
+  {
+    fthrow ( Exception, "The precomputation vector is zero...have you trained this classifier?" );
+  }
+
+  uint maxClassNo = 0;
+  for ( map<int, PrecomputedType>::const_iterator i = precomputedA.begin() ; i != precomputedA.end(); i++ )
+  {
+    uint classno = i->first;
+    maxClassNo = std::max ( maxClassNo, classno );
+    double beta;
+
+    if ( q != NULL ) {
+      map<int, double *>::const_iterator j = precomputedT.find ( classno );
+      double *T = j->second;
+      fmk->hik_kernel_sum_fast ( T, *q, xstar, beta );
+    } else {
+      const PrecomputedType & A = i->second;
+      map<int, PrecomputedType>::const_iterator j = precomputedB.find ( classno );
+      const PrecomputedType & B = j->second;
+
+      // fmk->hik_kernel_sum ( A, B, xstar, beta ); if A, B are of type Matrix
+      // Giving the transformation pf as an additional
+      // argument is necessary due to the following reason:
+      // FeatureMatrixT is sorted according to the original values, therefore,
+      // searching for upper and lower bounds ( findFirst... functions ) require original feature
+      // values as inputs. However, for calculation we need the transformed features values.
+
+      fmk->hik_kernel_sum ( A, B, xstar, beta, pf );
+    }
+
+    scores[ classno ] = beta;
+  }
+  scores.setDim ( maxClassNo + 1 );
+  
+  if ( precomputedA.size() > 1 ) {
+    // multi-class classification
+    return scores.maxElement();
+  } else {
+    // binary setting
+    // FIXME: not really flexible for every situation
+    scores[binaryLabelNegative] = -scores[binaryLabelPositive]; 
+    
+    return scores[ binaryLabelPositive ] <= 0.0 ? binaryLabelNegative : binaryLabelPositive;
+  }
+}
+
 void FMKGPHyperparameterOptimization::computePredictiveVarianceApproximateRough ( const NICE::SparseVector & x, NICE::Vector & predVariances ) const
 {
   double kSelf ( 0.0 );

+ 15 - 1
FMKGPHyperparameterOptimization.h

@@ -228,12 +228,26 @@ class FMKGPHyperparameterOptimization : NICE::Persistent
     /**
     * @brief classify an example 
     *
-    * @param x input example
+    * @param x input example (sparse vector)
     * @param scores scores for each class number
     *
     * @return class number achieving the best score
     */
     int classify ( const NICE::SparseVector & x, SparseVector & scores ) const;
+    
+    /**
+    * @brief classify an example that is given as non-sparse vector
+    * NOTE: whenever possible, you should sparse vectors to obtain significantly smaller computation times
+    * 
+    * @date 18-06-2013 (dd-mm-yyyy)
+    * @author Alexander Freytag
+    *
+    * @param x input example (non-sparse vector)
+    * @param scores scores for each class number
+    *
+    * @return class number achieving the best score
+    */
+    int classify ( const NICE::Vector & x, SparseVector & scores ) const;    
 
     /**
     * @brief compute predictive variance for a given test example using a rough approximation: k_{**} -  k_*^T (K+\sigma I)^{-1} k_* <= k_{**} - |k_*|^2 * 1 / \lambda_max(K + \sigma I), where we approximate |k_*|^2 by neglecting the mixed terms

+ 53 - 1
FastMinKernel.cpp

@@ -507,7 +507,6 @@ void FastMinKernel::hik_kernel_sum(const NICE::VVector & A, const NICE::VVector
 {
   // sparse version of hik_kernel_sum, no really significant changes,
   // we are just skipping zero elements
-  // for additional comments see the non-sparse version of hik_kernel_sum
   beta = 0.0;
   for (SparseVector::const_iterator i = xstar.begin(); i != xstar.end(); i++)
   {
@@ -559,6 +558,59 @@ void FastMinKernel::hik_kernel_sum(const NICE::VVector & A, const NICE::VVector
   }
 }
 
+void FastMinKernel::hik_kernel_sum(const NICE::VVector & A, const NICE::VVector & B, const NICE::Vector & xstar, double & beta, const ParameterizedFunction *pf) const
+{
+  beta = 0.0;
+  int dim ( 0 );
+  for (NICE::Vector::const_iterator i = xstar.begin(); i != xstar.end(); i++, dim++)
+  {
+  
+    double fval = *i;
+    
+    int nrZeroIndices = X_sorted.getNumberOfZeroElementsPerDimension(dim);
+    if ( nrZeroIndices == n ) {
+      // all features are zero and let us ignore it completely
+      continue;
+    }
+
+    int position;
+
+    //where is the example x^z_i located in
+    //the sorted array? -> perform binary search, runtime O(log(n))
+    // search using the original value
+    X_sorted.findFirstLargerInDimension(dim, fval, position);
+    position--;
+  
+    //NOTE again - pay attention! This is only valid if all entries are NOT negative! - if not, ask wether the current feature is greater than zero. If so, subtract the nrZeroIndices, if not do not
+    //sum_{l \in L_k} \alpha_l x^l_k
+    double firstPart(0.0);
+    //TODO in the "overnext" line there occurs the following error
+    // Invalid read of size 8
+    if (position >= 0) 
+      firstPart = (A[dim][position-nrZeroIndices]);
+    
+    // sum_{u \in U_k} alpha_u
+    
+    // sum_{u \in U_k} alpha_u
+    // => double secondPart( B(dim, n-1) - B(dim, position));
+    //TODO in the next line there occurs the following error
+    // Invalid read of size 8      
+    double secondPart( B[dim][n-1-nrZeroIndices]);
+    //TODO in the "overnext" line there occurs the following error
+    // Invalid read of size 8    
+    if (position >= 0) 
+      secondPart-= B[dim][position-nrZeroIndices];
+    
+    if ( pf != NULL )
+    {
+      fval = pf->f ( dim, fval );
+    }   
+    
+    // but apply using the transformed one
+    beta += firstPart + secondPart* fval;
+  }
+}
+
 void FastMinKernel::hik_kernel_sum_fast(const double *Tlookup, const Quantization & q, const NICE::Vector & xstar, double & beta) const
 {
   beta = 0.0;

+ 29 - 2
FastMinKernel.h

@@ -204,9 +204,25 @@ namespace NICE {
       */
       void hik_kernel_sum(const NICE::VVector & A, const NICE::VVector & B, const NICE::SparseVector & xstar, double & beta, const ParameterizedFunction *pf = NULL ) const;
       
+      /**
+      * @brief Computing k_{*}*alpha using the minimum kernel trick and exploiting sparsity of the feature vector given
+      * NOTE: Whenever possible, you should use sparse features to obtain significantly smaller computation times!
+      *
+      * @author Alexander Freytag
+      * @date 18-06-2013 (dd-mm-yyyy)
+      * @param A pre-computation matrix (VVector) (use the prepare method) 
+      * @param B pre-computation matrix (VVector)
+      * @param xstar new feature vector (non-sparse Vector)
+      * @param beta result of the scalar product
+      * @param pf optional feature transformation
+      */
+      void hik_kernel_sum(const NICE::VVector & A, const NICE::VVector & B, const NICE::Vector & xstar, double & beta, const ParameterizedFunction *pf = NULL ) const;      
+      
       /**
       * @brief compute beta = k_*^T * alpha by using a large lookup table created by hik_prepare_alpha_multiplications_fast
-      * @author Erik Rodner
+      * NOTE: Whenever possible, you should use sparse features to obtain significantly smaller computation times!
+      * @author Alexander Freytag
+      * @date 18-06-2013 (dd-mm-yyyy)
       *
       * @param Tlookup large lookup table calculated by hik_prepare_alpha_multiplications_fast
       * @param q Quantization object
@@ -214,11 +230,22 @@ namespace NICE {
       * @param beta result of the calculation
       */
       void hik_kernel_sum_fast(const double* Tlookup, const Quantization & q, const NICE::Vector & xstar, double & beta) const;
+      /**
+      * @brief compute beta = k_*^T * alpha by using a large lookup table created by hik_prepare_alpha_multiplications_fast
+      * NOTE: Whenever possible, you should use sparse features to obtain significantly smaller computation times!
+      * @author Alexander Frytag
+      *
+      * @param Tlookup large lookup table calculated by hik_prepare_alpha_multiplications_fast
+      * @param q Quantization object
+      * @param xstar feature vector (indirect k_*)
+      * @param beta result of the calculation
+      */      
+
       void hik_kernel_sum_fast(const double *Tlookup, const Quantization & q, const NICE::SparseVector & xstar, double & beta) const;
 
       /**
       * @brief compute lookup table for HIK calculation using quantized signals and prepare for K*alpha or k_*^T * alpha computations
-      * @author Erik Rodner
+      * @author Erik Rodner, Alexander Freytag
       *
       * @param alpha coefficient vector
       * @param A pre-calculation array computed by hik_prepare_alpha_multiplications

+ 1 - 1
FeatureMatrixT.h

@@ -81,7 +81,7 @@ template<class T> class FeatureMatrixT : NICE::Persistent
     /** 
     * @brief Recommended constructor
     * @author Alexander Freytag
-    * @date 07-12-2011 (dd-mm-yyyy)
+    * @date 07-12-2011 (dd-mm-yyyy) 
     */
     FeatureMatrixT(const std::vector<std::vector<T> > & _features, const int & _dim = -1);
     

+ 4 - 4
FeatureMatrixT.tcc

@@ -637,7 +637,7 @@ namespace NICE {
       {
         std::cerr << "set_features without features" << std::endl;
       }
-      
+            
       // resize our data structure      
       if (_dim >= 0) //did the user specified the number of dimensions?
         set_d(_dim);
@@ -664,8 +664,8 @@ namespace NICE {
             set_d(0);
           }          
         }
-      }    
-      
+      }  
+            
       // set number of examples n
       if (d>0)
       {
@@ -674,7 +674,7 @@ namespace NICE {
         else //we have examples x dimes (as usually done)   
           n = _features.size(); 
       }       
-      
+            
       // insert all values
       if (dimensionsOverExamples) //do we have dim x examples ?
       {

+ 38 - 0
GPHIKClassifier.cpp

@@ -116,6 +116,12 @@ void GPHIKClassifier::classify ( const SparseVector * example,  int & result, Sp
   this->classify( example, result, scores, tmpUncertainty );
 }
 
+void GPHIKClassifier::classify ( const NICE::Vector * example,  int & result, SparseVector & scores )
+{
+  double tmpUncertainty;
+  this->classify( example, result, scores, tmpUncertainty );
+}
+
 void GPHIKClassifier::classify ( const SparseVector * example,  int & result, SparseVector & scores, double & uncertainty )
 {
   scores.clear();
@@ -149,6 +155,38 @@ void GPHIKClassifier::classify ( const SparseVector * example,  int & result, Sp
   }    
 }
 
+void GPHIKClassifier::classify ( const NICE::Vector * example,  int & result, SparseVector & scores, double & uncertainty )
+{
+  scores.clear();
+  
+  int classno = gphyper->classify ( *example, scores );
+
+  if ( scores.size() == 0 ) {
+    fthrow(Exception, "Zero scores, something is likely to be wrong here: svec.size() = " << example->size() );
+  }
+  
+  result = scores.maxElement();
+   
+  if (uncertaintyPredictionForClassification)
+  {
+    if (varianceApproximation != NONE)
+    {
+      std::cerr << "ERROR: Uncertainty computation is currently not supported for NICE::Vector - use SparseVector instead" << std::endl;
+      uncertainty = std::numeric_limits<double>::max();
+    }  
+    else
+    {
+      //do nothing
+      uncertainty = std::numeric_limits<double>::max();
+    }
+  }
+  else
+  {
+    //do nothing
+    uncertainty = std::numeric_limits<double>::max();
+  }    
+}
+
 /** training process */
 void GPHIKClassifier::train ( const std::vector< NICE::SparseVector *> & examples, const NICE::Vector & labels )
 {

+ 23 - 0
GPHIKClassifier.h

@@ -89,6 +89,29 @@ class GPHIKClassifier
      * @param uncertainty (double*) predictive variance of the classification result, if computed
      */    
     void classify ( const NICE::SparseVector * example,  int & result, NICE::SparseVector & scores, double & uncertainty );
+    
+    /** 
+     * @brief classify a given example with the previously learnt model
+     * NOTE: whenever possible, you should the sparse version to obtain significantly smaller computation times* 
+     * @date 18-06-2013 (dd-mm-yyyy)
+     * @author Alexander Freytag
+     * @param example (non-sparse Vector) to be classified given in a non-sparse representation
+     * @param result (int) class number of most likely class
+     * @param scores (SparseVector) classification scores for known classes
+     */        
+    void classify ( const NICE::Vector * example,  int & result, NICE::SparseVector & scores );
+    
+    /** 
+     * @brief classify a given example with the previously learnt model
+     * NOTE: whenever possible, you should the sparse version to obtain significantly smaller computation times
+     * @date 18-06-2013 (dd-mm-yyyy)
+     * @author Alexander Freytag
+     * @param example (non-sparse Vector) to be classified given in a non-sparse representation
+     * @param result (int) class number of most likely class
+     * @param scores (SparseVector) classification scores for known classes
+     * @param uncertainty (double*) predictive variance of the classification result, if computed
+     */    
+    void classify ( const NICE::Vector * example,  int & result, NICE::SparseVector & scores, double & uncertainty );    
 
     /**
      * @brief train this classifier using a given set of examples and a given set of binary label vectors