classif.h File Reference

Include file for classification library. More...

#include <stdio.h>
#include <sys/types.h>
#include <signal.h>
#include <stdlib.h>
#include <math.h>
#include <time.h>
#include <string.h>
#include <gsl/gsl_rng.h>
#include <gsl/gsl_randist.h>
#include <gsl/gsl_statistics.h>
#include <misc.h>
Include dependency graph for classif.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Defines

#define _GNU_SOURCE
 GNU extensions.

Functions

void class_days_pc_clusters (int *days_class_cluster, double *pc_eof_days, double *eof_days_cluster, char *type, int neof, int ncluster, int ndays)
 Output a vector (dimension days) containing the closer (Euclidian distance) cluster number.
int generate_clusters (double *clusters, double *pc_eof_days, char *type, int nclassif, int neof, int ncluster, int ndays)
 Algorithm to generate clusters based on the Michelangeli et al (1995) methodology.
int best_clusters (double *best_clusters, double *pc_eof_days, char *type, int npart, int nclassif, int neof, int ncluster, int ndays)
 Algorithm to generate best clusters among many tries.
void mean_variance_dist_clusters (double *mean_dist, double *var_dist, double *pc, double *clusters, double *var_pc, double *var_pc_norm_all, int neof, int nclust, int ntime)
 Compute mean and variance of distances to clusters.
void dist_clusters_normctrl (double *dist_pc, double *pc, double *clusters, double *var_pc, double *var_pc_norm_all, double *mean_ctrl, double *var_ctrl, int neof, int nclust, int ntime)
 Compute distances to clusters normalized by control run mean and variance.

Detailed Description

Include file for classification library.

Definition in file classif.h.


Define Documentation

#define _GNU_SOURCE

GNU extensions.

Definition at line 55 of file classif.h.


Function Documentation

int best_clusters ( double *  best_clusters,
double *  pc_eof_days,
char *  type,
int  npart,
int  nclassif,
int  neof,
int  ncluster,
int  ndays 
)

Algorithm to generate best clusters among many tries.

Parameters:
[out] best_clusters Best clusters' positions.
[in] pc_eof_days Principal Components of EOF (daily data).
[in] type Type of distance used. Possible values: euclidian.
[in] npart Number of classification partitions to try.
[in] nclassif Maximum number of classifications to perform in the iterative algorithm.
[in] neof Number of EOFs.
[in] ncluster Number of clusters.
[in] ndays Number of days in the pc_eof_days vector.
Returns:
Minimum number of iterations needed.

Try to find best partition (clustering) which is closest to all the other partitions (which corresponds to the partition closest to the barycenter of partitions.

Definition at line 58 of file best_clusters.c.

References alloc_error(), and generate_clusters().

Referenced by main(), and wt_learning().

00058                                                                                                                                   {
00072   double min_meandistval; /* Minimum distance between a partition and all other partitions */
00073   double meandistval; /* Mean distance value between each corresponding clusters for the comparison of two partitions. */
00074   double maxdistval; /* Maximum distance over all clusters for the two partitions comparison. */
00075   double minval; /* Minimum distance to find a corresponding closest cluster in another partition. */
00076   double dist_bary; /* Distance summed over all EOFs between a cluster in one partition and other clusters in other partitions. */
00077   double val; /* Difference in positions between a cluster in one partition and other clusters in other partitions for a particular EOF. */
00078 
00079   double *tmpcluster = NULL; /* Temporary vector of clusters for one partition. */
00080   double *testclusters = NULL; /* Temporary vector of clusters for all partitions. */
00081   
00082   int min_cluster = -1; /* Cluster number used to find a corresponding cluster in another partition. */
00083   int min_partition = -1; /* Partition number used to find the partition which has the minimum distance to all other partitions. */
00084 
00085   int part; /* Loop counter for partitions */
00086   int part1; /* Loop counter for partitions inside loop */
00087   int part2; /* Loop counter for partitions inside loop */
00088   int clust; /* Loop counter for clusters */
00089   int clust1; /* Loop counter for clusters inside loop */
00090   int clust2; /* Loop counter for clusters inside loop */
00091   int eof; /* Loop counter for eofs */
00092 
00093   int niter; /* Number of iterations */
00094   int niter_min; /* Minimum number of iterations */
00095 
00096   (void) fprintf(stdout, "%s:: BEGIN: Find the best partition of clusters.\n", __FILE__);
00097 
00098   niter_min = 99999;
00099 
00100   /* Allocate memory */
00101   tmpcluster = (double *) calloc(neof*ncluster, sizeof(double));
00102   if (tmpcluster == NULL) alloc_error(__FILE__, __LINE__);
00103   testclusters = (double *) calloc(neof*ncluster*npart, sizeof(double));
00104   if (testclusters == NULL) alloc_error(__FILE__, __LINE__);
00105 
00106   /* Generate npart clusters (which will be used to find the best clustering). */
00107   (void) fprintf(stdout, "%s:: Generating %d partitions of clusters.\n", __FILE__, npart);
00108   for (part=0; part<npart; part++) {
00109 #if DEBUG >= 1
00110     (void) fprintf(stdout, "%s:: Generating %d/%d partition of clusters.\n", __FILE__, part+1, npart);
00111 #endif
00112     niter = generate_clusters(tmpcluster, pc_eof_days, type, nclassif, neof, ncluster, ndays);
00113     if (niter < niter_min) niter_min = niter;
00114     for (clust=0; clust<ncluster; clust++)
00115       for (eof=0; eof<neof; eof++)
00116         testclusters[part+eof*npart+clust*npart*neof] = tmpcluster[eof+clust*neof];
00117   }
00118 
00121   min_meandistval = 9999999999.9;
00122   min_partition = -1;
00123   /* Loop over all partition and compute distance between each other partition. */
00124   (void) fprintf(stdout, "%s:: Computing distance between each partitions of clusters.\n", __FILE__);
00125   for (part1=0; part1<npart; part1++) {
00126 #if DEBUG >= 1
00127     (void) fprintf(stdout, "%s:: Partition %d/%d.\n", __FILE__, part1+1, npart);
00128 #endif
00129     meandistval = 0.0;
00130     for (part2=0; part2<npart; part2++) {
00131 
00132       /* Don't compute for the same partition number. */
00133       if (part1 != part2) {
00134 
00135         maxdistval = -9999999999.9;
00136         
00137         for (clust1=0; clust1<ncluster; clust1++) {
00138           
00139           /* Find closest cluster to current one (in terms of distance summed over all EOF). */
00140           minval = 9999999999.9;
00141           min_cluster = -1;
00142           for (clust2=0; clust2<ncluster; clust2++) {
00143 
00144             if ( !strcmp(type, "euclidian") ) {
00145               /* Sum distances over all EOF. */
00146               dist_bary = 0.0;
00147               for (eof=0; eof<neof; eof++) {
00148                 val = testclusters[part2+eof*npart+clust1*npart*neof] - testclusters[part1+eof*npart+clust2*npart*neof];
00149                 dist_bary += (val * val);
00150               }
00151               
00152               dist_bary = sqrt(dist_bary);
00153             }
00154             else {
00155               (void) fprintf(stderr, "best_clusters: ABORT: Unknown distance type=%s!!\n", type);
00156               (void) abort();
00157             }
00158             
00159             /* Check for minimum distance. We want to find the corresponding closest cluster in another partition. */
00160             if (dist_bary < minval) {
00161               minval = dist_bary;
00162               min_cluster = clust2;
00163             }
00164           }
00165 
00166           if (min_cluster == -1) {
00167             (void) fprintf(stderr, "best_clusters: ABORT: Error in algorithm. Cannot find best cluster!\n");
00168             (void) abort();
00169           }
00170           
00171           /* Save the maximum distance over all clusters for the two partitions comparison. */
00172           if (minval > maxdistval)
00173             maxdistval = minval;
00174         }
00175         /* Sum the maximum distance of the clusters between each corresponding one over all the partitions.
00176            We want to compute the mean afterward. */
00177         meandistval += maxdistval;
00178       }
00179     }
00180     /* Compute the mean of the distances between each corresponding clusters for the comparison of two partitions. */
00181     meandistval = meandistval / (double) (npart-1);
00182     /* We want to keep the partition which has the minimum distance to all other partitions. */
00183     if (meandistval < min_meandistval) {
00184       min_meandistval = meandistval;
00185       min_partition = part1;
00186     }
00187   }
00188 
00189   if (min_partition == -1) {
00190     /* Failing algorithm */
00191     (void) fprintf(stderr, "best_clusters: ABORT: Error in algorithm. Cannot find best partition!\n");
00192     (void) abort();
00193   }
00194 
00195   /* Save data for the best selected partition of clusters. */
00196   (void) fprintf(stdout, "%s:: Save best partition of clusters.\n", __FILE__);
00197   for (clust=0; clust<ncluster; clust++)
00198     for (eof=0; eof<neof; eof++)
00199       best_clusters[eof+clust*neof] = testclusters[min_partition+eof*npart+clust*npart*neof];  
00200 
00201   /* Free memory. */
00202   (void) free(tmpcluster);
00203   (void) free(testclusters);
00204 
00205   (void) fprintf(stdout, "%s:: END: Find the best partition of clusters. Partition %d selected.\n", __FILE__, min_partition);
00206 
00207   return niter_min;
00208 }

void class_days_pc_clusters ( int *  days_class_cluster,
double *  pc_eof_days,
double *  eof_days_cluster,
char *  type,
int  neof,
int  ncluster,
int  ndays 
)

Output a vector (dimension days) containing the closer (Euclidian distance) cluster number.

The distance is computed as the distance between a day's principal components for each EOF and the cluster principal components for each EOF. To evaluate the closest one, all the square of the distances for all EOFs are summed for each cluster, before the square root is applied.

Parameters:
[out] days_class_cluster Cluster number associated for each day.
[in] pc_eof_days Principal Components of EOF (daily data).
[in] eof_days_cluster Clusters' centroid positions for each eof.
[in] type Type of distance used. Possible values: euclidian.
[in] neof Number of EOFs.
[in] ncluster Number of clusters.
[in] ndays Number of days in the pc_eof_days vector.

Definition at line 60 of file class_days_pc_clusters.c.

Referenced by generate_clusters(), wt_downscaling(), and wt_learning().

00061                                                           {
00072   double dist_min; /* Minimum distance found between a given day PC (summed over all EOF) and each cluster centroid. */
00073   int clust_dist_min; /* Cluster number which has the minimum distance dist_min */
00074   double dist_sum; /* Sum of distances (partial computation) over all EOFs */
00075   double val; /* Distance between a given day PC (for a particular EOF) and one cluster centroid. */
00076   double dist_clust; /* Distance (full computation of dist_sum). */
00077   
00078   int day; /* Loop counter for days */
00079   int clust; /* Loop counter for cluster */
00080   int eof; /* Loop counter for eofs */
00081 
00082   if ( !strcmp(type, "euclidian") ) {
00083     /* Euclidian distance type */
00084 
00085     /* Parse each day */
00086     for (day=0; day<ndays; day++) {
00087 
00088       /* Initialize */
00089       dist_min = 9999999999.0;
00090       clust_dist_min = 999;
00091 
00092 #if DEBUG >= 7
00093       (void) fprintf(stderr, "day=%d\n", day);
00094 #endif
00095 
00096       /* Parse each cluster */
00097       for (clust=0; clust<ncluster; clust++) {
00098 
00099 #if DEBUG >= 7
00100         (void) fprintf(stderr, "clust=%d\n", clust);
00101 #endif
00102 
00103         dist_sum = 0.0;
00104         /* Sum all distances (over EOF) between the PC of the day and the PC of the cluster centroid for each EOF respectively */
00105         for (eof=0; eof<neof; eof++) {
00106           val = pc_eof_days[day+eof*ndays] - eof_days_cluster[eof+clust*neof];
00107 #if DEBUG >= 9
00108           printf("%d %d %lf %lf\n",clust,eof,pc_eof_days[day+eof*ndays],eof_days_cluster[eof+clust*neof]);
00109 #endif
00110           /* Euclidian distance: square */
00111           dist_sum += (val * val);
00112         }
00113         /* Euclidian distance: square root of squares */
00114         dist_clust = sqrt(dist_sum);
00115 
00116 #if DEBUG >= 7
00117         (void) fprintf(stderr, "dist_clust=%lf\n", dist_clust);
00118 #endif
00119 
00120         /* Is it a cluster which has less distance as the minimum found yet ? */
00121         if (dist_clust < dist_min) {
00122           /* Save cluster number */
00123           clust_dist_min = clust;
00124           dist_min = dist_clust;
00125         }
00126       }
00127       if (clust_dist_min == 999) {
00128         /* Failing algorithm */
00129         (void) fprintf(stderr, "%s: ABORT: Impossible: no cluster was selected!! Problem in algorithm...\n", __FILE__);
00130         (void) abort();
00131       }
00132       /* Assign cluster with minimum distance to all EOFs for this day */
00133       days_class_cluster[day] = clust_dist_min;
00134 #if DEBUG >= 9
00135       (void) fprintf(stderr, "%s: day %d cluster %d\n", __FILE__, day, clust_dist_min);
00136 #endif
00137     }
00138   }
00139   else {
00140     /* Unknown distance type */
00141     (void) fprintf(stderr, "%s: ABORT: Unknown distance type=%s!!\n", __FILE__, type);
00142     (void) abort();
00143   }
00144 }

void dist_clusters_normctrl ( double *  dist_pc,
double *  pc,
double *  clusters,
double *  var_pc,
double *  var_pc_norm_all,
double *  mean_ctrl,
double *  var_ctrl,
int  neof,
int  nclust,
int  ntime 
)

Compute distances to clusters normalized by control run mean and variance.

Parameters:
[out] dist_pc Distances (normalized by control run mean and variance) of normalized EOF-projected large-scale field to clusters
[in] pc EOF-projected large-scale field
[in] clusters Cluster centroids for each EOF in EOF-projected space of the large-scale field
[in] var_pc Variance of EOF-projected large-scale field of the learning period, for each EOF separately.
[in] var_pc_norm_all Norm of the variance of the first EOF of the EOF-projected large-scale field of the control run.
[in] mean_ctrl Mean of the distances to clusters for the control run, for each cluster separately.
[in] var_ctrl Variance of the distances to clusters for the control run, for each cluster separately.
[in] neof EOF dimension
[in] nclust Clusters dimension
[in] ntime Time dimension

Definition at line 59 of file dist_clusters_normctrl.c.

Referenced by wt_downscaling(), and wt_learning().

00061                                                         {
00075   double sum;
00076   double val;
00077   
00078   int eof;
00079   int nt;
00080   int clust;
00081 
00082   for (clust=0; clust<nclust; clust++) {
00083     for (nt=0; nt<ntime; nt++) {
00084       sum = 0.0;
00085       for (eof=0; eof<neof; eof++) {
00086         val = (pc[nt+eof*ntime] / sqrt(var_pc_norm_all[eof])) - (clusters[eof+clust*neof] / sqrt(var_pc[eof]));
00087         sum += (val * val);
00088         //        if (nt == (ntime-5)) {
00089         //          printf("dist_val %d %d %d %lf %lf %lf %lf %lf\n",clust,nt,eof,val,pc[nt+eof*ntime],var_pc_norm_all[eof],clusters[eof+clust*neof],var_pc[eof]);
00090         //        }
00091       }
00092       dist_pc[nt+clust*ntime] = ( (sqrt(sum) - mean_ctrl[clust]) / sqrt(var_ctrl[clust]) );
00093       //      if (nt == (ntime -5))
00094       //        printf("dist_pc %d %d %lf %lf %lf %lf %lf\n",clust,nt,dist_pc[nt+clust*ntime],sqrt(sum),mean_ctrl[clust],var_ctrl[clust]);
00095       //    }
00096     }
00097   }
00098 }

int generate_clusters ( double *  clusters,
double *  pc_eof_days,
char *  type,
int  nclassif,
int  neof,
int  ncluster,
int  ndays 
)

Algorithm to generate clusters based on the Michelangeli et al (1995) methodology.

Parameters:
[out] clusters Clusters' positions.
[in] pc_eof_days Principal Components of EOF (daily data).
[in] type Type of distance used. Possible values: euclidian.
[in] nclassif Maximum number of classifications to perform in the iterative algorithm.
[in] neof Number of EOFs.
[in] ncluster Number of clusters.
[in] ndays Number of days in the pc_eof_days vector.
Returns:
Number of iterations.

Generate ncluster random days

Warning: we are using time() as the seed. Don't run this subroutine twice with the same time. If you do you will get the exact same sequence.

Main algorithm

Try to find the maximum distance (PC-space) between the new cluster center and the previous value

Definition at line 58 of file generate_clusters.c.

References alloc_error(), and class_days_pc_clusters().

Referenced by best_clusters(), and main().

00059                                                      {
00072   unsigned long int *random_num = NULL; /* Vector of random numbers for random choice of initial points. */
00073   int eof; /* Loop counter for EOF. */
00074   int clust; /* Loop counter for clusters. */
00075   int day; /* Loop counter for days. */
00076   int classif; /* Loop counter for classifications. */
00077   int ndays_cluster; /* Number of days in the current cluster. */
00078 
00079   const gsl_rng_type *T; /* For random number generation type. */
00080   gsl_rng *rng; /* For random number generation. */
00081 
00082   double cluster_bary; /* Cluster barycentre. */
00083   double ndiff_cluster_bary; /* Distance between current cluster barycenter and previous value in the iteration. */
00084   double mean_days; /* Mean of the days (PC-space) for a cluster. */
00085   double *eof_days_cluster = NULL; /* Vector of clusters' barycenter positions (PC-space). */
00086   int *days_class_cluster = NULL; /* Vector of classification of days into each cluster. */
00087 
00088   static unsigned long int seed = 0;
00089 
00090   (void) fprintf(stdout, "%s:: BEGIN: Find clusters among data points.\n", __FILE__);
00091 
00092   /***********************************/
00094   /***********************************/
00095 
00096   (void) fprintf(stdout, "%s:: Choosing %d random points.\n", __FILE__, ncluster);
00097 
00098   /* Initialize random number generator */
00099   T = gsl_rng_default;
00100   rng = gsl_rng_alloc(T);
00103   if (seed == 0) seed = time(NULL);
00104   (void) gsl_rng_set(rng, seed++);
00105   
00106   /* Generate ncluster random days and initialize cluster PC array */
00107   random_num = (unsigned long int *) calloc(ncluster, sizeof(unsigned long int));
00108   if (random_num == NULL) alloc_error(__FILE__, __LINE__);
00109   for (clust=0; clust<ncluster; clust++)
00110     random_num[clust] = gsl_rng_uniform_int(rng, ndays);  
00111   
00112   /* Free random number generator */
00113   (void) gsl_rng_free(rng);
00114   
00115   /********************/
00117   /********************/
00118 
00119   /* Allocate memory */
00120   eof_days_cluster = (double *) calloc(neof*ncluster, sizeof(double));
00121   if (eof_days_cluster == NULL) alloc_error(__FILE__, __LINE__);
00122   days_class_cluster = (int *) calloc(ndays, sizeof(int));
00123   if (days_class_cluster == NULL) alloc_error(__FILE__, __LINE__);
00124   
00125   /* Initialize cluster PC array randomly */
00126   (void) fprintf(stdout, "%s:: Initializing cluster array.\n", __FILE__);
00127   for (eof=0; eof<neof; eof++)
00128     for (clust=0; clust<ncluster; clust++) {
00129       eof_days_cluster[eof+clust*neof] = pc_eof_days[random_num[clust]+eof*ndays];
00130 
00131 #if DEBUG >= 7
00132       (void) fprintf(stderr, "eof=%d cluster=%d eof_days_cluster=%lf\n", eof, clust, eof_days_cluster[eof+clust*neof]);
00133 #endif
00134     }
00135 
00136   (void) free(random_num);
00137 
00138   /* Iterate by performing up to nclassif classifications. Stop if same cluster center positions in two consecutive iterations. */
00139   cluster_bary = -9999999999.9;
00140   ndiff_cluster_bary = -9999999999.9;
00141   classif = 0;
00142 
00143   (void) fprintf(stdout, "%s:: Iterate up to %d classifications or when classification is stable.\n", __FILE__, nclassif);
00144 
00145   while (ndiff_cluster_bary != 0.0 && classif < nclassif) {
00146 
00147 #if DEBUG >= 7
00148     (void) fprintf(stderr, "classif=%d cluster_bary=%lf\n", classif, cluster_bary);
00149 #endif
00150 
00151     /* Classify each day (pc_eof_days) in the current clusters (eof_days_cluster) = days_class_cluster */
00152     (void) class_days_pc_clusters(days_class_cluster, pc_eof_days, eof_days_cluster, type, neof, ncluster, ndays);
00153 
00154     /* For each cluster, perform a mean of all points falling in that cluster.
00155        Compare to the current clusters by calculating the 'coordinates' (PC-space) of the 'new' cluster center. */
00156     cluster_bary = -9999999999.9;
00157 
00158     /* Loop over clusters and EOFs */
00159     for (clust=0; clust<ncluster; clust++) {
00160       for (eof=0; eof<neof; eof++) {
00161         mean_days = 0.0;
00162         ndays_cluster = 0;
00163 
00164         /* Compute the mean (PC-space) of all the days of the current cluster */
00165         for (day=0; day<ndays; day++)
00166           if (days_class_cluster[day] == clust) {
00167             mean_days += pc_eof_days[day+eof*ndays];
00168             ndays_cluster++;
00169           }
00170 
00171         if (ndays_cluster > 0) {
00172 
00173           mean_days = mean_days / (double) ndays_cluster;
00174 
00177 #if DEBUG >= 7
00178           (void) fprintf(stderr, "eof=%d cluster=%d diff_cluster_bary=%lf mean_days=%lf eof_days_cluster=%lf\n",
00179                          eof, clust, fabs(mean_days - eof_days_cluster[eof+clust*neof]), mean_days, eof_days_cluster[eof+clust*neof]);
00180 #endif
00181 
00182           /* Compute the difference between the new cluster center position and the previous one */
00183           ndiff_cluster_bary = fabs(mean_days - eof_days_cluster[eof+clust*neof]);
00184 
00185           /* If this new cluster center position is further away than the other EOF's ones, choose this new cluster center */
00186           if ( ndiff_cluster_bary > cluster_bary )
00187             cluster_bary = mean_days;
00188 
00189           /* Store the new cluster center value */
00190           clusters[eof+clust*neof] = mean_days;
00191 
00192 #if DEBUG >= 7
00193           if (classif == 0 || cluster_bary == 0.0)
00194             (void) fprintf(stderr, "eof=%d cluster=%d mean_pc_days=%lf ndays_cluster=%d cluster_bary=%lf\n",
00195                            eof, clust, mean_days, ndays_cluster, cluster_bary);
00196 #endif
00197         }
00198         else
00199           clusters[eof+clust*neof] = 0;
00200       }
00201     }
00202 
00203     classif++;
00204 
00205     /* Update the cluster center matrix with the new values */
00206     if (cluster_bary != 0.0 && classif < nclassif)
00207       for (clust=0; clust<ncluster; clust++)
00208         for (eof=0; eof<neof; eof++)
00209           eof_days_cluster[eof+clust*neof] = clusters[eof+clust*neof];
00210   }
00211 
00212   /* Free memory */
00213   (void) free(eof_days_cluster);
00214   (void) free(days_class_cluster);
00215 
00216   (void) fprintf(stdout, "%s:: END: Find clusters among data points. %d iterations needed.\n", __FILE__, classif);
00217 
00218   return classif;
00219 }

void mean_variance_dist_clusters ( double *  mean_dist,
double *  var_dist,
double *  pc,
double *  clusters,
double *  var_pc,
double *  var_pc_norm_all,
int  neof,
int  nclust,
int  ntime 
)

Compute mean and variance of distances to clusters.

Parameters:
[out] mean_dist Mean of distances to clusters.
[out] var_dist Variance of distances to clusters.
[in] pc EOF-projected large-scale field.
[in] clusters Cluster centroids for each EOF in EOF-projected space of the large-scale field.
[in] var_pc Norm of the variance of the first EOF of the EOF-projected large-scale field of the learning period.
[in] var_pc_norm_all Norm of the variance of the first EOF of the EOF-projected large-scale field of the control run.
[in] neof EOF dimension
[in] nclust Clusters dimension
[in] ntime Time dimension

Definition at line 58 of file mean_variance_dist_clusters.c.

References alloc_error().

Referenced by main(), and wt_downscaling().

00059                                                                                       {
00072   double *dist_pc = NULL; /* Sum of the normalized distances to clusters over all EOFs */
00073 
00074   double sum = 0.0; /* Sum of normalized distances to clusters */
00075   double val; /* Normalized distance to cluster */
00076   
00077   int eof; /* EOF loop counter */
00078   int nt; /* Time loop counter */
00079   int clust; /* Cluster loop counter */
00080 
00081   /* Allocate memory */
00082   dist_pc = (double *) malloc(ntime * sizeof(double));
00083   if (dist_pc == NULL) alloc_error(__FILE__, __LINE__);
00084 
00085   /* Loop over all clusters */
00086   for (clust=0; clust<nclust; clust++) {
00087     /* Loop over time */
00088     for (nt=0; nt<ntime; nt++) {
00089       /* Initialize the sum */
00090       sum = 0.0;
00091       /* Loop over all EOFs */
00092       for (eof=0; eof<neof; eof++) {
00093         /* Calculate normalized distance to clusters */
00094         val = (pc[nt+eof*ntime] / sqrt(var_pc_norm_all[eof])) - (clusters[eof+clust*neof] / sqrt(var_pc[eof]));
00095         /* Calculate squared sum */
00096         sum += (val * val);
00097       }
00098       /* Save for this timestep the square root of the sum of the squares of the distance to clusters */
00099       dist_pc[nt] = sqrt(sum);
00100     }
00101     /* Calculate mean over time */
00102     mean_dist[clust] = gsl_stats_mean(dist_pc, 1, ntime);
00103     /* Calculate variance over time */
00104     var_dist[clust] = gsl_stats_variance(dist_pc, 1, ntime);
00105 
00106     //    printf("ctrl dist pre... %d %lf %lf %lf %lf %lf\n",clust,sqrt(sum),pc[ntime-1+clust*ntime],sqrt(var_pc_norm_all[clust]),clusters[clust+clust*neof],sqrt(var_pc[clust]));
00107     //    printf("ctrl dist.... %d %lf %lf %lf %lf\n",clust,dist_pc[ntime-1],sum,mean_dist[clust],sqrt(var_dist[clust]));
00108   }
00109 
00110   /* Free memory */
00111   (void) free(dist_pc);
00112 }


Generated on 12 May 2016 for DSCLIM by  doxygen 1.6.1