openModeller  Version 1.4.0
environmental_distance.cpp
Go to the documentation of this file.
00001 //
00002 // Generic environmental distance algorithm
00003 //
00004 // Description: Generic algorithm based on distances.
00005 //
00006 // Authors:      Danilo J. S. Bellini <danilo.estagio@gmail.com>
00007 //               Renato De Giovanni <renato@cria.org.br>
00008 // Copyright:   See COPYING file that comes with this distribution
00009 // Date:        2006-09-18
00010 //
00011 
00012 #include "environmental_distance.hh"
00013 #include <openmodeller/Exceptions.hh>
00014 #include <openmodeller/ScaleNormalizer.hh>
00015 
00016 //
00017 // METADATA
00018 //
00019 #define DATA_MIN 0.0
00020 #define DATA_MAX 1.0
00021 #define NUM_PARAM 3
00022 
00023 #define PARDISTTYPE "DistanceType"
00024 #define PARPOINTQNT "NearestPoints"
00025 #define PARDIST     "MaxDistance"
00026 #define PARDISTMIN 0.0
00027 #define PARDISTMAX 1.0
00028 
00029 #define LOG_SQRT_PI 0.5723649429247000870717135 /* log (sqrt (pi)) */
00030 #define I_SQRT_PI   0.5641895835477562869480795 /* 1 / sqrt (pi) */
00031 #define BIGX        20.0  /* max value to represent exp (x) */
00032 #define ex(x)       (((x) < -BIGX) ? 0.0 : exp (x))
00033 #define Z_MAX       6.0   /* maximum meaningful z value */
00034 
00035 static AlgParamMetadata parameters[NUM_PARAM] = { // Parameters
00036    { // 1st parameter
00037       PARDISTTYPE, // Id
00038       "Metric",    // Name
00039       Integer,     // Type
00040       "Metric used to calculate distances: " // Overview
00041          "1=Euclidean, "
00042          "2=Mahalanobis, "
00043          "3=Manhattan/Gower, "
00044          "4=Chebyshev.",
00045       "Metric used to calculate distances: " // Description
00046          "1=Euclidean, "
00047          "2=Mahalanobis, "
00048          "3=Manhattan/Gower, "
00049          "4=Chebyshev",
00050       1,                                               // Not zero if the parameter has lower limit
00051       FIRST_DISTANCE_TYPE,                             // Parameter's lower limit
00052       1,                                               // Not zero if the parameter has upper limit
00053       FIRST_DISTANCE_TYPE + AMOUNT_DISTANCE_TYPES - 1, // Parameter's upper limit
00054       "1"         // Parameter's typical (default) value
00055    },
00056    { // 2nd parameter
00057       PARPOINTQNT,            // Id
00058       "Nearest \'n\' points", // Name
00059       Integer,                // Type
00060       "Nearest \'n\' points whose mean value will be the reference when calculating environmental distances.", // Overview
00061       "Nearest \'n\' points whose mean value will be the reference when calculating environmental distances. When set to 1, distances will be measured to the closest point, which is the same behavior of the previously existing minimum distance algorithm. When set to 0, distances will be measured to the average of all presence points, which is the same behavior of the previously existing distance to average algorithm. Intermediate values between 1 and the total number of presence points are now accepted.", // Description
00062       1,          // Not zero if the parameter has lower limit
00063       0,          // Parameter's lower limit
00064       1,          // Not zero if the parameter has upper limit
00065       1000,       // Parameter's upper limit
00066       "1"         // Parameter's typical (default) value
00067    },
00068    { // 3rd parameter
00069       PARDIST,            // Id
00070       "Maximum distance", // Name
00071       Real,               // Type
00072       "Maximum distance to the reference in the environmental space.", // Overview
00073       "Maximum distance to the reference in the environmental space, above which the conditions will be considered unsuitable for presence. Since 1 corresponds to the biggest possible distance between any two points in the environment space, setting the maximum distance to this value means that all points in the environmental space will have an associated probability. The probability of presence for points that fall within the range of the maximum distance is inversely proportional to the distance to the reference point (linear decay). The only exception is when the maximum distance is 1 and the metric is Mahalanobis, which will produce probabilities following the chi-square distribution.", // Description
00074       1,          // Not zero if the parameter has lower limit
00075       PARDISTMIN, // Parameter's lower limit
00076       1,          // Not zero if the parameter has upper limit
00077       PARDISTMAX, // Parameter's upper limit
00078       "0.1"       // Parameter's typical (default) value
00079    },
00080 };
00081 
00082 static AlgMetadata metadata = { // General metadata
00083   "ENVDIST",                    // Id
00084   "Environmental Distance",     // Name
00085   "0.5",                        // Version
00086   "Generic algorithm based on environmental dissimilarity metrics.", // Overview
00087   "Generic algorithm based on environmental dissimilarity metrics. When used with the Gower metric and maximum distance 1, this algorithm should produce the same result of the algorithm known as DOMAIN.", // Description
00088   "Mauro E. S. Munoz, Renato De Giovanni, Danilo J. S. Bellini",    // Algorithm author
00089   "Carpenter G, Gillison AN, Winter J (1993) DOMAIN: A flexible modeling procedure for mapping potential distributions of animals and plants. Biodiversity and Conservation 2: 667-680. Farber O & Kadmon R 2003. Assessment of alternative approaches for bioclimatic modeling with special emphasis on the Mahalanobis distance. Ecological Modelling 160: 115–130.", // Bibliography
00090   "Danilo J. S. Bellini, Renato De Giovanni", // Code author
00091   "danilo . estagio [at] gmail . com, renato [at] cria . org . br", // Code author's contact
00092   0,                    // Does not accept categorical data
00093   0,                    // Does not need (pseudo)absence points
00094   NUM_PARAM, parameters // Algorithm's parameters
00095 };
00096 
00097 
00098 //
00099 // LINK TO OM
00100 //
00101 
00102 // Needed code to link this algorithm to oM
00103 OM_ALG_DLL_EXPORT
00104 AlgorithmImpl *algorithmFactory(){
00105    return new EnvironmentalDistance(); // Create an instance of the algorithm's class
00106 }
00107 OM_ALG_DLL_EXPORT
00108 AlgMetadata const *algorithmMetadata(){
00109    return &metadata;
00110 }
00111 
00112 
00113 //
00114 // ALGORITHM CONSTRUCTOR/DESTRUCTOR
00115 //
00116 
00117 // Constructor for the algorithm class
00118 EnvironmentalDistance::EnvironmentalDistance() : AlgorithmImpl(&metadata){
00119    _cov_matrix = _cov_matrix_inv = NULL;
00120    _normalizerPtr = new ScaleNormalizer( DATA_MIN, DATA_MAX, true );
00121 }
00122 
00123 // Destructor for the algorithm class
00124 EnvironmentalDistance::~EnvironmentalDistance()
00125 {
00126    switch(_par_dist_type){
00127       case MahalanobisDistance:
00128          if(_cov_matrix!=NULL){
00129             delete _cov_matrix;
00130             delete _cov_matrix_inv;
00131          }
00132          break;
00133       //case ManhattanDistance:
00134       //case ChebyshevDistance:
00135       //case EuclideanDistance:
00136       //default:
00137    }
00138 }
00139 
00140 
00141 //
00142 // ALGORITHM GENERIC METHODS (virtual AlgorithmImpl methods)
00143 //
00144 
00145 // Initialize the algorithm
00146 int EnvironmentalDistance::initialize(){
00147 
00148    // Test the parameters' data types
00149    if(!getParameter(PARDIST,&_par_dist)){
00150       Log::instance()->error("Parameter '" PARDIST "' was not passed.\n");
00151       return 0;
00152    }
00153    if(!getParameter(PARDISTTYPE,&_par_dist_type)){
00154       Log::instance()->error("Parameter '" PARDISTTYPE "' was not passed.\n");
00155       return 0;
00156    }
00157    if(!getParameter(PARPOINTQNT,&_par_point_qnt)){
00158       Log::instance()->error("Parameter '" PARPOINTQNT "' was not passed.\n");
00159       return 0;
00160    }
00161 
00162    // Impose limits to the parameters, if somehow the user don't obey
00163    if     (_par_dist>PARDISTMAX) _par_dist = PARDISTMAX;
00164    else if(_par_dist<PARDISTMIN) _par_dist = PARDISTMIN;
00165 
00166    // Check distance type
00167    switch(_par_dist_type){
00168       case EuclideanDistance:
00169          Log::instance()->debug("Using Euclidean distance\n");
00170          break;
00171       case MahalanobisDistance:
00172          Log::instance()->debug("Using Mahalanobis distance\n");
00173          break;
00174       case ManhattanDistance:
00175          Log::instance()->debug("Using Manhattan distance\n");
00176          break;
00177       case ChebyshevDistance:
00178          Log::instance()->debug("Using Chebyshev distance\n");
00179          break;
00180       default:
00181          Log::instance()->error("Parameter '" PARDISTTYPE "' wasn't set properly. It should be an integer between 1 and 4.\n");
00182          return 0;
00183    }
00184 
00185    //_samp->environmentallyUnique(); // Debug
00186 
00187    // Initialize some common-use attributes
00188    _layer_count    = _samp->numIndependent();
00189    _presence_count = _samp->numPresence();
00190 
00191    // Load all environmental data of presence points
00192    if(_presence_count == 0){
00193       Log::instance()->error("There is no presence point.\n");
00194       return 0;
00195    }
00196 
00197    // Check number of nearest points parameter
00198    if (_par_point_qnt > _presence_count)
00199       Log::instance()->warn("Parameter '" PARPOINTQNT "' is greater than the number of presence points\n");
00200    else if (_par_point_qnt < 0) _par_point_qnt = 0;
00201 
00202    OccurrencesPtr presences = _samp->getPresences();
00203    for(int i = 0 ; i < _presence_count ; i++)
00204       _presence_points.push_back((*presences)[i]->environment());
00205    // Calcs the mean of all presence points
00206    _average_point = _presence_points[0]; // There is at least one presence point
00207    for(int i = 1 ; i < _presence_count ; i++)
00208       _average_point += _presence_points[i];
00209    _average_point /= _presence_count;
00210 
00211    _use_chisq = false;
00212 
00213    // Allow using "Distance" method and normalize _par_dist
00214    if(!_init_distance_type()){
00215       Log::instance()->error("Could not determine maximum distance in the environmental space.\n");
00216       return 0;
00217    }
00218 
00219    _done = true;       // Needed for not-iterative algorithms
00220    return 1; // There was no problem in initialization
00221 }
00222 
00223 // Returns the occurrence probability
00224 Scalar EnvironmentalDistance::getValue(const Sample& x) const{
00225    Scalar dist;
00226 
00227    //
00228    // Distance to average
00229    //
00230    if((_par_point_qnt >= _presence_count) || (_par_point_qnt <= 0))
00231       dist = _distance(x, _average_point);
00232 
00233    //
00234    // Minimum distance (not really needed)
00235    //
00236    else if(_par_point_qnt == 1){
00237       Scalar distIterator;
00238       dist = -1;
00239       for(int i = 0 ; i < _presence_count ; i++){ // Iterate for each presence point
00240          distIterator = _distance(x, _presence_points[i]);
00241          if((distIterator < dist || dist < 0)){
00242             dist = distIterator;
00243          }
00244       }
00245 
00246    //
00247    // Mean of _par_point_qnt nearest points
00248    //
00249    }else{
00250       Scalar distIterator, distTmp;
00251       int indexIterator, indexTmp;
00252       Sample nearMean;
00253       std::vector<int> nearestIndex(_par_point_qnt);
00254       std::vector<Scalar> nPdist(_par_point_qnt);
00255       //Log::instance()->debug("Starting with distances:\n"); // Debug
00256       for(int i = 0 ; i < _par_point_qnt ; i++){ // We know that _par_point_qnt < _presence_count
00257          nPdist[i] = _distance(x, _presence_points[i]);
00258          //x.dump(); // debug
00259          //_presence_points[i].dump(); // debug
00260          nearestIndex[i] = i;
00261          //Log::instance()->debug("   dist[%d]=%.8g\n", nearestIndex[i], nPdist[i]); // Debug
00262       }
00263 
00264       //Log::instance()->debug("Nearest points:\n"); // Debug
00265       for(int i = _par_point_qnt ; i < _presence_count ; i++){ // This loop finds the nearest points
00266          distIterator = _distance(x, _presence_points[i]);
00267          indexIterator = i;
00268          //Log::instance()->debug("dist[%d] = %.8g:\n", indexIterator, distIterator); // Debug
00269          for(int j = 0 ; j < _par_point_qnt ; j++){ // Trade pointIterator with the first smaller point
00270             if(nPdist[j] > distIterator){
00271                distTmp = distIterator;
00272                indexTmp = indexIterator;
00273                distIterator = nPdist[j];
00274                indexIterator = nearestIndex[j];
00275                nPdist[j] = distTmp;
00276                nearestIndex[j] = indexTmp;
00277                j = -1; // Re-start the for loop for the new value
00278             }
00279          }
00280          //for(int j = 0 ; j < _par_point_qnt ; j++) // Debug
00281          //   Log::instance()->debug("   dist[%d]=%.8g\n", nearestIndex[j], nPdist[j]);
00282       }
00283 
00284       // Now we have the nearest points. Let's get its mean:
00285       nearMean = _presence_points[nearestIndex[0]]; // There is at least one point
00286       for(int i = 1 ; i < _par_point_qnt ; i++)
00287          nearMean += _presence_points[nearestIndex[i]];
00288       nearMean /= _par_point_qnt;
00289 
00290       dist = _distance(x, nearMean);
00291    }
00292 
00293    //Log::instance()->debug("distance=%.8g\n\n", dist); // Debug
00294    //Log::instance()->debug("max dist=%.8g\n\n", _par_dist); // Debug
00295 
00296    // Now finishes the algorithm calculating the probability
00297    if(dist < 0) // There isn't any occurrence
00298       return 0.0;
00299    else if(_use_chisq) // Only for Mahalanobis distance when maxdist == 1
00300       return _pochisq(dist, _layer_count-1 );
00301    else if(dist > _par_dist) // Point is too faraway from nearest point
00302       return 0.0;
00303    else
00304       return 1.0 - (dist / _par_dist);
00305 }
00306 
00307 // Initialize _cov_matrix and _cov_matrix_inv
00308 void EnvironmentalDistance::_calc_covariance_matrix(){
00309    if(_cov_matrix!=NULL){ // Garbage collector
00310       delete _cov_matrix;
00311       delete _cov_matrix_inv;
00312    }
00313    _cov_matrix = new Matrix(_layer_count,_layer_count); // Alloc memory for new matrices
00314    _cov_matrix_inv = new Matrix(_layer_count,_layer_count);
00315 
00316    if(_presence_count > 1){
00317       // Calcs the cross-covariance for each place in the matrix
00318       for(int i = 0 ; i < _layer_count ; i++)
00319          for(int j = i ; j < _layer_count ; j++){
00320             (*_cov_matrix)(i,j) = 0;
00321             for(int k = 0 ; k < _presence_count ; k++)
00322                (*_cov_matrix)(i,j) += (_presence_points[k][i]-_average_point[i]) *
00323                                     (_presence_points[k][j]-_average_point[j]);
00324             (*_cov_matrix)(i,j) /= _presence_count;
00325             (*_cov_matrix)(j,i) = (*_cov_matrix)(i,j);
00326          }
00327    }
00328    else{
00329       // Use an identity matrix if there's only one point
00330       for(int i = 0 ; i < _layer_count ; i++)
00331          for(int j = i ; j < _layer_count ; j++){
00332             (*_cov_matrix)(i,j) = (i == j) ? 1 : 0;
00333          }
00334    }
00335    //Log::instance()->debug("Cov matrix:\n");
00336    //std::cout << (*_cov_matrix); // Debug
00337 
00338    try{
00339       (*_cov_matrix_inv) = !(*_cov_matrix);
00340    }
00341    catch ( std::exception& e ) {
00342       string msg = e.what();
00343       msg.append( "\nExperiment has no solution using Mahalanobis distance.\n" );
00344       Log::instance()->error( "%s", msg.c_str() );
00345       throw AlgorithmException( msg.c_str() );
00346    }
00347    //Log::instance()->debug("Inv cov matrix:\n");
00348    //std::cout << (*_cov_matrix_inv); // Debug
00349 }
00350 
00351 // Calcs the distance between x and y using _par_dist_type
00352 inline Scalar EnvironmentalDistance::_distance(const Sample& x, const Sample& y) const{
00353    Scalar dist=0;
00354    switch(_par_dist_type){
00355 
00356       //
00357       // Mahalanobis Distance
00358       //
00359       case MahalanobisDistance:{
00360          Matrix lineMatrix(1,_layer_count);
00361          // Make lineMatrix = x - y
00362          for(int i=0; i<_layer_count; i++)
00363             lineMatrix(0,i) = x[i] - y[i];
00364          // Definition of Mahalanobis distance
00365          dist = (lineMatrix * (*_cov_matrix_inv) * (~lineMatrix))(0,0); // Operator () of a 1x1 matrix
00366          if(!_use_chisq)
00367             dist = sqrt(dist);
00368          //Log::instance()->info("\nDISTANCE: %g\n",dist); // Debug
00369       }break;
00370 
00371       //
00372       // Manhattan Distance
00373       //
00374       case ManhattanDistance:{
00375          Scalar tmp;
00376          for(int k=0;k<_layer_count;k++){
00377             tmp = x[k] - y[k];
00378             // We don't need this because we did that in _par_dist normalization:
00379             // tmp /= DATA_MAX - DATA_MIN; // range(k)
00380             if(tmp < 0)
00381                dist -= tmp;
00382             else
00383                dist += tmp;
00384          }
00385          dist /= _layer_count;
00386       }break;
00387 
00388       //
00389       // ChebyshevDistance
00390       //
00391       case ChebyshevDistance:
00392       {
00393          Scalar tmp;
00394          for(int i=0; i<_layer_count; i++){
00395             tmp = x[i] - y[i];
00396             if(tmp < 0)
00397               tmp = -tmp;
00398             if(tmp > dist)
00399               dist = tmp;
00400          }
00401       }break;
00402 
00403       //
00404       // Euclidean Distance
00405       //
00406       case EuclideanDistance:
00407       default:{
00408          Sample delta = x;
00409          delta -= y;
00410          dist = delta.norm(); // The usual norm of the "delta" vector have
00411       }break;                 // Euclidean distance for a n-dimensional space
00412    }
00413 
00414    return dist;
00415 }
00416 
00417 // Initialize the data structures of a distance type and normalize _par_dist
00418 bool EnvironmentalDistance::_init_distance_type(){
00419    Scalar distMax;
00420 
00421    // Calcs the maximum distance (distMax)
00422    switch(_par_dist_type){
00423 
00424       case MahalanobisDistance:{
00425          _calc_covariance_matrix(); // Initialize _cov_matrix_inv
00426          if(_par_dist < 1.0) {
00427             Log::instance()->info("Using normalized maximum distance\n"); // Debug
00428             Scalar distIterator;
00429             Sample x,y;
00430             x.resize(_layer_count);
00431             y.resize(_layer_count);
00432             distMax = 0.0;
00433             bool foundDist = false;
00434             // Distance between oposite edges (vertex) of the hypercube
00435             for(int i=0; i<(1<<(_layer_count-1)); i++){ // for(i FROM 0 TO 2 "power" (_layer_count - 1))
00436                for(int k=0; k<_layer_count; k++) // This is the same loop used to create
00437                   if((i & (1<<k)) != 0){         // binary numbers, but with "max" and
00438                      x[k] = DATA_MAX;          // "min" instead of "1" and "0"
00439                      y[k] = DATA_MIN; // y is the oposite of x
00440                   }else{
00441                      x[k] = DATA_MIN;
00442                      y[k] = DATA_MAX;
00443                   }
00444                distIterator = _distance(x,y);
00445                if(distIterator > distMax){
00446                   distMax = distIterator;
00447                   foundDist = true;
00448                }
00449             }
00450             if(!foundDist)
00451                return false;
00452          }
00453          else {
00454             // In this case, chi-square probabilities will be used, 
00455             // so no need to find the max distance
00456             Log::instance()->info("Using chi-square probabilities\n");
00457             _use_chisq = true;
00458             return true;
00459          }
00460       }break;
00461 
00462       case ManhattanDistance:
00463          // Manhattan maximum value is always DATA_MAX - DATA_MIN, even for n
00464          // dimensions, so there's no real initialization for Manhattan
00465       case ChebyshevDistance:
00466          // Chebyshev and Manhattan distances has the same maximum value
00467          distMax = DATA_MAX - DATA_MIN;
00468          break;
00469 
00470       case EuclideanDistance:
00471       default:
00472          //   For a 1 dimensional world, we impose a size limit of L
00473          //   For a 2D world, we have a square, with size limit of L*sqrt(2)
00474          // since each side is a 1 dimensional world with max = 1
00475          //   In a 3D world, size limit is L*sqrt(3) because we have a cube
00476          //   For a n-dimensional world, we have the maximum dist equals
00477          // L*sqrt(n), because sqrt((L*sqrt(n-1))^2 + L^2) = L*sqrt(n)
00478          // (indution using Pythagoras).
00479          distMax = sqrt((double)_layer_count) * (DATA_MAX - DATA_MIN);
00480    }
00481 
00482    // Normalize _par_dist for its limits
00483    _par_dist = (_par_dist - PARDISTMIN)/(PARDISTMAX - PARDISTMIN); // Now _par_dist is in [0,1]
00484    _par_dist *= distMax; // Now _par_dist is in [0,distMax]
00485    return true;
00486 }
00487 
00488 Scalar EnvironmentalDistance::_poz(Scalar z) const {
00489    Scalar  y, x, w;
00490    if (z == 0.0)
00491       x = 0.0;
00492    else {
00493       y = 0.5 * fabs (z);
00494       if (y >= (Z_MAX * 0.5))
00495          x = 1.0;
00496       else if (y < 1.0) {
00497          w = y*y;
00498          x = ((((((((0.000124818987 * w
00499                     -0.001075204047) * w +0.005198775019) * w
00500                     -0.019198292004) * w +0.059054035642) * w
00501                     -0.151968751364) * w +0.319152932694) * w
00502                     -0.531923007300) * w +0.797884560593) * y * 2.0;
00503       }
00504       else {
00505          y -= 2.0;
00506          x = (((((((((((((-0.000045255659 * y
00507                           +0.000152529290) * y -0.000019538132) * y
00508                           -0.000676904986) * y +0.001390604284) * y
00509                           -0.000794620820) * y -0.002034254874) * y
00510                           +0.006549791214) * y -0.010557625006) * y
00511                           +0.011630447319) * y -0.009279453341) * y
00512                           +0.005353579108) * y -0.002141268741) * y
00513                           +0.000535310849) * y +0.999936657524;
00514        }
00515     }
00516   return (z > 0.0 ? ((x + 1.0) * 0.5) : ((1.0 - x) * 0.5));
00517 }
00518 
00519 Scalar EnvironmentalDistance::_pochisq(Scalar x, int df) const {
00520    Scalar  a, y = 0, s;
00521    Scalar  e, c, z;
00522    int     even;     /* true if df is an even number */
00523    if (x <= 0.0 || df < 1)
00524       return (0.0);
00525    a = 0.5 * x;
00526    even = (2*(df/2)) == df;
00527    if (df > 1)
00528       y = ex (-a);
00529    s = (even ? y : (2.0 * _poz(-sqrt (x))));
00530    if (df > 2) {
00531       x = 0.5 * (df - 1.0);
00532       z = (even ? 1.0 : 0.5);
00533       if (a > BIGX) {
00534          e = (even ? 0.0 : LOG_SQRT_PI);
00535          c = log (a);
00536          while (z <= x) {
00537             e = log (z) + e;
00538             s += ex (c*z-a-e);
00539             z += 1.0;
00540          }
00541          return (s);
00542       }
00543       else {
00544          e = (even ? 1.0 : (I_SQRT_PI / sqrt(a)));
00545          c = 0.0;
00546          while (z <= x) {
00547             e = e * (a / z);
00548             c = c + e;
00549             z += 1.0;
00550          }
00551          return (c * y + s);
00552       }
00553    }
00554    else
00555       return (s);
00556 }
00557 
00558 // Alg serializer
00559 void EnvironmentalDistance::_getConfiguration(ConfigurationPtr& config) const {
00560    if (!_done ) return;
00561 
00562    ConfigurationPtr model_config( new ConfigurationImpl("EnvironmentalDistance") );
00563    config->addSubsection(model_config);
00564    model_config->addNameValue("MaxDistance",_par_dist);
00565    ConfigurationPtr envpoints_config(new ConfigurationImpl("EnvironmentalReferences"));
00566    model_config->addSubsection(envpoints_config);
00567    for(unsigned int i=0; i<_presence_points.size(); i++) {
00568       ConfigurationPtr point_config(new ConfigurationImpl("Reference"));
00569       envpoints_config->addSubsection(point_config);
00570       point_config->addNameValue("Value", _presence_points[i]);
00571    }
00572    model_config->addNameValue("Average",_average_point);
00573 }
00574 
00575 // Alg deserializer
00576 void EnvironmentalDistance::_setConfiguration(const ConstConfigurationPtr& config) {
00577    ConstConfigurationPtr model_config = config->getSubsection("EnvironmentalDistance",false);
00578 
00579    if (!model_config) return;
00580 
00581    // Metric
00582    if(!getParameter(PARDISTTYPE,&_par_dist_type)){
00583       Log::instance()->error("Parameter '" PARDISTTYPE "' was not found in serialized model.\n");
00584       return;
00585    }
00586    // "n" closest points
00587    if(!getParameter(PARPOINTQNT,&_par_point_qnt)){
00588       Log::instance()->error("Parameter '" PARPOINTQNT "' was not found in serialized model.\n");
00589       return;
00590    }
00591    Scalar max_distance;
00592    if(!getParameter(PARDIST,&max_distance)){
00593       Log::instance()->error("Parameter '" PARDIST "' was not found in serialized model.\n");
00594       return;
00595    }
00596    _use_chisq = false;
00597    if (_par_dist_type == MahalanobisDistance && max_distance == 1.0) {
00598       _use_chisq = true;
00599    }
00600    // Maximum distance
00601    _par_dist = model_config->getAttributeAsDouble("MaxDistance", 0.0);
00602    // Environmental points
00603    ConstConfigurationPtr envpoints_config = model_config->getSubsection("EnvironmentalReferences",false);
00604    Configuration::subsection_list subs = envpoints_config->getAllSubsections();
00605    Configuration::subsection_list::iterator begin = subs.begin();
00606    Configuration::subsection_list::iterator end = subs.end();
00607    for (; begin != end; ++begin) {
00608       if ((*begin)->getName() != "Reference") continue;
00609       Sample point = (*begin)->getAttributeAsSample("Value");
00610       _presence_points.push_back(point);
00611    }
00612    // Average
00613    _average_point = model_config->getAttributeAsSample("Average");
00614 
00615    _layer_count = (int)_average_point.size();
00616    _presence_count = (int)_presence_points.size();
00617 
00618    if ( _par_dist_type == MahalanobisDistance ) {
00619       _calc_covariance_matrix(); // Initialize _cov_matrix_inv
00620    }
00621 
00622    _done = true;
00623 }