Main Page   Modules   Data Structures   File List   Data Fields   Globals   Related Pages  

gnn_prototype.c

Go to the documentation of this file.
00001 /***************************************************************************
00002  *  @file gnn_prototype.c
00003  *  @brief Prototype Implementation File.
00004  *
00005  *  @date   : 17-09-03 20:52
00006  *  @author : Pedro Ortega C. <peortega@dcc.uchile.cl>
00007  *  Copyright  2003  Pedro Ortega C.
00008  ****************************************************************************/
00009 /*
00010  *  This program is free software; you can redistribute it and/or modify
00011  *  it under the terms of the GNU General Public License as published by
00012  *  the Free Software Foundation; either version 2 of the License, or
00013  *  (at your option) any later version.
00014  *
00015  *  This program is distributed in the hope that it will be useful,
00016  *  but WITHOUT ANY WARRANTY; without even the implied warranty of
00017  *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00018  *  GNU Library General Public License for more details.
00019  *
00020  *  You should have received a copy of the GNU General Public License
00021  *  along with this program; if not, write to the Free Software
00022  *  Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
00023  */
00024 
00025 
00026 
00027 /**
00028  * @defgroup gnn_prototype_doc gnn_prototype : Prototypes.
00029  * @ingroup gnn_atomic_doc
00030  *
00031  * The \ref gnn_prototype datatype implements a set of prototypes in
00032  * \f$\mathbb{R}^n\f$, which are compared against the input vectors.
00033  * The function computed by this node is:
00034  *
00035  * \f[ y_j = \sqrt{ \sum_{i=1}^n (w^j_i - x_i)^2 } \f]
00036  *
00037  * where \f$w^j\f$ is the \f$j\f$-th prototype. In other words, the outputs
00038  * \f$y_j\f$ measure the euclidian distance to each one of the prototypes.
00039  *
00040  * The function's gradients are:
00041  *
00042  * \f[ \frac{\partial E}{\partial x_i}
00043  *     = - \frac{1}{2} \left (
00044  *           \sum_{i=1}^n (w_i^j - x_i)^2
00045  *                     \right )^{-\frac{1}{2}} 2 (w_i^j - x_i)
00046  *         \frac{\partial E}{\partial y}
00047  *     = - \frac{w_i^k - x_i}{y_j} \frac{\partial E}{\partial y}
00048  * \f]
00049  * and
00050  * \f[ \frac{\partial E}{\partial w^j_i}
00051  *     = \frac{1}{2} \left (
00052  *           \sum_{i=1}^n (w_i^j - x_i)^2
00053  *                     \right )^{-\frac{1}{2}} 2 (w_i^j - x_i)
00054  *       \frac{\partial E}{\partial y}
00055  *     = \frac{w_i^k - x_i}{y_j} \frac{\partial E}{\partial y}
00056  * \f]
00057  * The prototypes are very effective in cojunction with a kernel function,
00058  * like \ref gnn_gaussian, thus forming a Kernel Model.
00059  */
00060 
00061 
00062 
00063 /******************************************/
00064 /* Include Files                          */
00065 /******************************************/
00066 
00067 #include "gnn_prototype.h"
00068 #include "gnn_utilities.h"
00069 #include <gsl/gsl_blas.h>
00070 #include <math.h>
00071 
00072 
00073 
00074 /******************************************/
00075 /* Static Declaration                     */
00076 /******************************************/
00077 
00078 static int
00079 gnn_prototype_f (gnn_node *node,
00080                  const gsl_vector *x,
00081                  const gsl_vector *w,
00082                  gsl_vector *y);
00083 
00084 static int
00085 gnn_prototype_dx (gnn_node *node,
00086                   const gsl_vector *x,
00087                   const gsl_vector *w,
00088                   const gsl_vector *dy,
00089                   gsl_vector *dx);
00090 
00091 static int
00092 gnn_prototype_dw (gnn_node *node,
00093                   const gsl_vector *x,
00094                   const gsl_vector *w,
00095                   const gsl_vector *dy,
00096                   gsl_vector *dw);
00097 
00098 static void
00099 gnn_prototype_destroy (gnn_node *node);
00100 
00101 
00102 
00103 /******************************************/
00104 /* Static Implementation                  */
00105 /******************************************/
00106 
00107 /**
00108  * @brief Computes the output.
00109  * @ingroup gnn_prototype_doc
00110  *
00111  * This functions evaluates the distance to the prototypes.
00112  *
00113  * @param node  A pointer to a \ref gnn_prototype node.
00114  * @param x The input vector \f$x\f$.
00115  * @param w The parameter vector \f$w\f$.
00116  * @param y The output buffer vector \f$y\f$.
00117  * @return Returns 0 on success.
00118  */
00119 static int
00120 gnn_prototype_f (gnn_node *node,
00121                  const gsl_vector *x,
00122                  const gsl_vector *w,
00123                  gsl_vector *y)
00124 {
00125     size_t j;
00126     size_t m;
00127     gnn_prototype *pnode;
00128     gsl_vector_view wj;
00129 
00130     /* get the sizes */
00131     m = gnn_node_output_get_size (node);
00132 
00133     /* get gnn_prototype view */
00134     pnode = (gnn_prototype *) node;
00135 
00136     /* evaluate */
00137     for (j=0; j<m; ++j)
00138     {
00139         double yj;
00140 
00141         /* compute euclidean distance to prototype */
00142 /*
00143         wj = gsl_matrix_row (pnode->D, j);
00144         gsl_vector_memcpy (pnode->row, &(wj.vector));
00145         gsl_vector_sub (pnode->row, x);
00146         yj = gsl_blas_dnrm2 (pnode->row);
00147 */
00148 
00149         wj = gsl_matrix_row (pnode->D, j);
00150         yj = gnn_vector_euclidian_dist (&(wj.vector), x);
00151 
00152         /* set distance */
00153         gsl_vector_set (pnode->y, j, yj);
00154     }
00155 
00156     /* sum result from buffer */
00157     gsl_vector_memcpy (y, pnode->y);
00158 
00159     return 0;
00160 }
00161 
00162 /**
00163  * @brief Computes \f$ \frac{\partial E}{\partial x} \f$.
00164  * @ingroup gnn_prototype_doc
00165  *
00166  * This functions computes the gradient of the distance to the prototypes
00167  * with respect to the inputs. That is,
00168  * \f[ \frac{\partial E}{\partial x_i}
00169  *     = - \sum_{j=1}^m \frac{w^j_i - x_i}{y_j} \frac{\partial E}{\partial y_j}
00170  * \f]
00171  *
00172  * @param node  A pointer to a \ref gnn_prototype node.
00173  * @param x  The input vector \f$x\f$.
00174  * @param w  The parameter vector \f$w\f$.
00175  * @param dy The vector \f$\frac{\partial E}{\partial y}\f$.
00176  * @param dx The output buffer vector \f$\frac{\partial E}{\partial x}\f$.
00177  * @return Returns 0 on success.
00178  */
00179 static int
00180 gnn_prototype_dx (gnn_node *node,
00181                   const gsl_vector *x,
00182                   const gsl_vector *w,
00183                   const gsl_vector *dy,
00184                   gsl_vector *dx)
00185 {
00186     size_t i;
00187     size_t n;
00188     gnn_prototype *pnode;
00189     gsl_vector_view wi;
00190 
00191     /* get the size */
00192     n = gnn_node_input_get_size (node);
00193 
00194     /* get gnn_prototype view */
00195     pnode = (gnn_prototype *) node;
00196 
00197     /* evaluate */
00198     for (i=0; i<n; ++i)
00199     {
00200         double xi;
00201         double dxi;
00202 
00203         /* get i-th component */
00204         xi  = gsl_vector_get (x, i);
00205         //dxi = gsl_vector_get (dx, i);
00206 
00207         /* compute dyj/dxi */
00208         gsl_vector_set_all (pnode->col, xi);
00209         wi = gsl_matrix_column (pnode->D, i);
00210         gsl_vector_sub (pnode->col, &(wi.vector));
00211         gsl_vector_div (pnode->col, pnode->y);
00212 
00213         /* compute dE/dxi */
00214         gsl_vector_mul (pnode->col, dy);
00215         //dxi += gnn_vector_sum_elements (pnode->col);
00216         dxi = gnn_vector_sum_elements (pnode->col);
00217 
00218         /* set dxi */
00219         gsl_vector_set (dx, i, dxi);
00220     }
00221 
00222     return 0;
00223 }
00224 
00225 /**
00226  * @brief Computes \f$ \frac{\partial E}{\partial w} \f$.
00227  * @ingroup gnn_prototype_doc
00228  *
00229  * This functions computes the gradient of the distance to the prototypes
00230  * with respect to the prototypes' position. That is,
00231  * \f[ \frac{\partial E}{\partial w^j_i}
00232  *     = \frac{w^j_i - x_i}{y_j} \frac{\partial E}{\partial y_j}
00233  * \f]
00234  *
00235  * @param node  A pointer to a \ref gnn_prototype node.
00236  * @param x  The input vector \f$x\f$.
00237  * @param w  The parameter vector \f$w\f$.
00238  * @param dy The vector \f$\frac{\partial E}{\partial y}\f$.
00239  * @param dw The output buffer vector \f$\frac{\partial E}{\partial w}\f$.
00240  * @return Returns 0 on success.
00241  */
00242 static int
00243 gnn_prototype_dw (gnn_node *node,
00244                   const gsl_vector *x,
00245                   const gsl_vector *w,
00246                   const gsl_vector *dy,
00247                   gsl_vector *dw)
00248 {
00249     size_t j;
00250     size_t m;
00251     gnn_prototype *pnode;
00252     gsl_vector_view wj;
00253     gsl_vector_view dwj;
00254 
00255     /* get the sizes */
00256     m = gnn_node_output_get_size (node);
00257 
00258     /* get gnn_prototype view */
00259     pnode = (gnn_prototype *) node;
00260 
00261     /* evaluate (compute row-wise) */
00262     for (j=0; j<m; ++j)
00263     {
00264         double yj;
00265         double dyj;
00266 
00267         /* get output y and dy */
00268         yj  = gsl_vector_get (pnode->y, j);
00269         dyj = gsl_vector_get (dy, j);
00270 
00271         /* compute dE/dwj */
00272         wj = gsl_matrix_row (pnode->D, j);
00273         gsl_vector_memcpy (pnode->row, &(wj.vector));
00274         gsl_vector_sub (pnode->row, x);
00275         gsl_vector_scale (pnode->row, dyj / yj);
00276 
00277         /* set dE/dwj */
00278         dwj = gsl_matrix_row (pnode->dD, j);
00279         gsl_vector_add (&(dwj.vector), pnode->row);
00280     }
00281 
00282     return 0;
00283 }
00284 
00285 /**
00286  * @brief The \ref gnn_prototype destroy function.
00287  * @ingroup gnn_prototype_doc
00288  *
00289  * This functions is the \ref gnn_prototype destroy function.
00290  *
00291  * @param  node A node to an \ref gnn_prototype.
00292  */
00293 static void
00294 gnn_prototype_destroy (gnn_node *node)
00295 {
00296     gnn_prototype *pnode;
00297 
00298     assert (node != NULL);
00299     
00300     pnode = (gnn_prototype *) node;
00301     
00302     if (pnode->y != NULL)
00303         gsl_vector_free (pnode->y);
00304     if (pnode->row != NULL)
00305         gsl_vector_free (pnode->row);
00306     if (pnode->col != NULL)
00307         gsl_vector_free (pnode->col);
00308 
00309     return;
00310 }
00311 
00312 
00313 
00314 /******************************************/
00315 /* Public Interface                       */
00316 /******************************************/
00317 
00318 /**
00319  * @brief Creates a Prototype node.
00320  * @ingroup gnn_prototype_doc
00321  *
00322  * This function creates a node of the \ref gnn_prototype type. Each one
00323  * of its vectors is drawn from a uniform distribution between
00324  * [min; max).
00325  *
00326  * @param input_size  The input size \f$n\f$.
00327  * @param output_size The output size \f$m\f$, and the number of prototypes.
00328  * @param min         The minimum value of a component.
00329  * @param max         The maximum value of a component.
00330  * @return A pointer to a new \ref gnn_prototype node.
00331  */
00332 gnn_node *
00333 gnn_prototype_new (size_t input_size, size_t output_size,
00334                                               double min, double max)
00335 {
00336     int status;
00337     size_t i, j;
00338     gsl_rng  *rng;
00339     gnn_node *node;
00340     gnn_prototype *pnode;
00341     gsl_vector *w;
00342     gsl_vector *dw;
00343     gsl_vector_int *f;
00344 
00345     /* check if sizes are positive */
00346     if (input_size < 1)
00347     {
00348         GSL_ERROR_VAL ("input size should be strictly positive",
00349                        GSL_EINVAL, NULL);
00350     }
00351     if (output_size < 1)
00352     {
00353         GSL_ERROR_VAL ("output size should be strictly positive",
00354                        GSL_EINVAL, NULL);
00355     }
00356 
00357     /* check interval */
00358     if (min >= max)
00359     {
00360         GSL_ERROR_VAL ("min should be smaller than max",
00361                        GSL_EINVAL, NULL);
00362     }
00363 
00364     /* allocate the node */
00365     pnode = (gnn_prototype *) malloc (sizeof (gnn_prototype));
00366     if (pnode == NULL)
00367     {
00368         GSL_ERROR_VAL ("couldn't allocate memory for gnn_prototype",
00369                        GSL_ENOMEM, NULL);
00370     }
00371     
00372     /* get view as gnn_node */
00373     node = (gnn_node *) pnode;
00374 
00375     /* Initialize the node */
00376     status = gnn_node_init (node,
00377                             "gnn_prototype",
00378                             gnn_prototype_f,
00379                             gnn_prototype_dx,
00380                             gnn_prototype_dw,
00381                             gnn_prototype_destroy);
00382     if (status)
00383     {
00384         free (node);
00385         GSL_ERROR_VAL ("could not initialize gnn_prototype node",
00386                        GSL_EINVAL, NULL);
00387     }
00388 
00389     /* set sizes */
00390     status = gnn_node_set_sizes (node,  input_size,
00391                                         output_size, input_size * output_size);
00392     if (status)
00393     {
00394         gnn_node_destroy (node);
00395         GSL_ERROR_VAL ("could not set sizes for gnn_prototype node",
00396                        GSL_EFAILED, NULL);
00397     }
00398 
00399     /* Allocate memory for buffers */
00400     pnode->y   = gsl_vector_alloc (output_size);
00401     pnode->row = gsl_vector_alloc (input_size);
00402     pnode->col = gsl_vector_alloc (output_size);
00403     if (pnode->y == NULL || pnode->row == NULL || pnode->col == NULL)
00404     {
00405         gnn_node_destroy (node);
00406         GSL_ERROR_VAL ("could not allocate buffers for gnn_prototype",
00407                        GSL_ENOMEM, NULL);
00408     }
00409 
00410     /* Get parameter views */
00411     w  = gnn_node_local_get_w (node);
00412     dw = gnn_node_local_get_dw (node);
00413     f  = gnn_node_local_get_f (node);
00414     
00415     pnode->D_view  = gsl_matrix_view_vector (w, output_size, input_size);
00416     pnode->dD_view = gsl_matrix_view_vector (dw, output_size, input_size);
00417     pnode->Df_view = gsl_matrix_int_view_vector (f, output_size, input_size);
00418     
00419     pnode->D  = &(pnode->D_view.matrix);
00420     pnode->dD = &(pnode->dD_view.matrix);
00421     pnode->Df = &(pnode->Df_view.matrix);
00422 
00423     /* initialize prototypes */
00424     rng = gnn_get_rng ();
00425     for (j=0; j<output_size; ++j)
00426         for (i=0; i<input_size; ++i)
00427         {
00428             double w;
00429             w = gsl_rng_uniform (rng);
00430             w = (max - min) * w + min;
00431             gsl_matrix_set (pnode->D, j, i, w);
00432         }
00433 
00434     /* update parameters */
00435     gnn_node_local_update (node);
00436 
00437     return node;
00438 }
00439 
00440 /**
00441  * @brief Creates a standard Prototype node.
00442  * @ingroup gnn_prototype_doc
00443  *
00444  * This function creates a node of the \ref gnn_prototype type. Each one
00445  * of its vectors is drawn from a uniform distribution between
00446  * [-1; 1).
00447  *
00448  * @param input_size  The input size \f$n\f$.
00449  * @param output_size The output size \f$m\f$, and the number of prototypes.
00450  * @return A pointer to a new \ref gnn_prototype node.
00451  */
00452 gnn_node *
00453 gnn_prototype_standard_new (size_t input_size, size_t output_size)
00454 {
00455     return gnn_prototype_new (input_size, output_size, -1.0, 1.0);
00456 }
00457 
00458 /**
00459  * @brief Freeze a prototype.
00460  * @ingroup gnn_prototype_doc
00461  *
00462  * This function freezes the j-th prototype.
00463  *
00464  * @param node A pointer to a \ref gnn_prototype node.
00465  * @param j    The index of the prototype to be frozen.
00466  * @return Returns 0 if succeeded.
00467  */
00468 int
00469 gnn_prototype_vector_freeze (gnn_node *node, size_t j)
00470 {
00471     gnn_prototype *pnode;
00472     gsl_vector_int_view v;
00473 
00474     assert (node != NULL);
00475 
00476     /* get view */
00477     pnode = (gnn_prototype *) node;
00478 
00479     /* check for boundaries */
00480     if (j < 0 || j >= pnode->Df->size1)
00481     {
00482         GSL_ERROR ("index out of bounds", GSL_EINVAL);
00483     }
00484 
00485     /* get freeze flags */
00486     gnn_node_local_get_f (node);
00487 
00488     /* get vector view */
00489     v = gsl_matrix_int_row (pnode->Df, j);
00490 
00491     /* freeze vector */
00492     gsl_vector_int_set_all (&(v.vector), 1);
00493 
00494     /* update parameters */
00495     gnn_node_local_update (node);
00496 
00497     return 0;
00498 }
00499 
00500 /**
00501  * @brief Unfreeze a prototype.
00502  * @ingroup gnn_prototype_doc
00503  *
00504  * This function unfreezes the j-th prototype.
00505  *
00506  * @param node A pointer to a \ref gnn_prototype node.
00507  * @param j    The index of the prototype to be unfrozen.
00508  * @return Returns 0 if succeeded.
00509  */
00510 int
00511 gnn_prototype_vector_unfreeze (gnn_node *node, size_t j)
00512 {
00513     gnn_prototype *pnode;
00514     gsl_vector_int_view v;
00515 
00516     assert (node != NULL);
00517 
00518     /* get view */
00519     pnode = (gnn_prototype *) node;
00520 
00521     /* check for boundaries */
00522     if (j < 0 || j >= pnode->Df->size1)
00523     {
00524         GSL_ERROR ("index out of bounds", GSL_EINVAL);
00525     }
00526 
00527     /* get freeze flags */
00528     gnn_node_local_get_f (node);
00529 
00530     /* get vector view */
00531     v = gsl_matrix_int_row (pnode->Df, j);
00532 
00533     /* freeze vector */
00534     gsl_vector_int_set_zero (&(v.vector));
00535 
00536     /* update parameters */
00537     gnn_node_local_update (node);
00538 
00539     return 0;
00540 }
00541 
00542 /**
00543  * @brief Sets a prototype.
00544  * @ingroup gnn_prototype_doc
00545  *
00546  * This function sets a new value for the j-th prototype.
00547  *
00548  * @param node A pointer to a \ref gnn_prototype node.
00549  * @param j    The index of the prototype to be frozen.
00550  * @param v    A new with the new values for the prototype.
00551  * @return Returns 0 if succeeded.
00552  */
00553 int
00554 gnn_prototype_vector_set (gnn_node *node, size_t j, const gsl_vector *v)
00555 {
00556     gsl_vector_view u;
00557     gnn_prototype *pnode;
00558 
00559     assert (node != NULL);
00560 
00561     /* get view */
00562     pnode = (gnn_prototype *) node;
00563 
00564     /* check for boundaries */
00565     if (j < 0 || j >= pnode->D->size1)
00566     {
00567         GSL_ERROR ("index out of bounds", GSL_EINVAL);
00568     }
00569 
00570     /* check size */
00571     if (gnn_node_input_get_size (node) != v->size)
00572     {
00573         GSL_ERROR ("vector insn't of the correct size", GSL_EINVAL);
00574     }
00575 
00576     /* get parameters */
00577     gnn_node_local_get_w (node);
00578 
00579     /* get vector view */
00580     u = gsl_matrix_row (pnode->D, j);
00581 
00582     /* copy vector */
00583     gsl_vector_memcpy (&(u.vector), v);
00584 
00585     /* update parameters */
00586     gnn_node_local_update (node);
00587 
00588     return 0;
00589 }
00590 
00591 /**
00592  * @brief Gets a prototype.
00593  * @ingroup gnn_prototype_doc
00594  *
00595  * This function gets the values of the j-th prototype and stores them
00596  * into the given vector.
00597  *
00598  * @param node A pointer to a \ref gnn_prototype node.
00599  * @param j    The index of the prototype to be frozen.
00600  * @param v    A vector of the correct size, where the prototype's values
00601  *             should be copied in.
00602  * @return Returns 0 if succeeded.
00603  */
00604 int
00605 gnn_prototype_vector_get (gnn_node *node, size_t j, gsl_vector *v)
00606 {
00607     gsl_vector_view u;
00608     gnn_prototype *pnode;
00609 
00610     assert (node != NULL);
00611 
00612     /* get view */
00613     pnode = (gnn_prototype *) node;
00614 
00615     /* check for boundaries */
00616     if (j < 0 || j >= pnode->D->size1)
00617     {
00618         GSL_ERROR ("index out of bounds", GSL_EINVAL);
00619     }
00620 
00621     /* check size */
00622     if (gnn_node_input_get_size (node) != v->size)
00623     {
00624         GSL_ERROR ("vector insn't of the correct size", GSL_EINVAL);
00625     }
00626 
00627     /* get parameters */
00628     gnn_node_local_get_w (node);
00629 
00630     /* get vector view */
00631     u = gsl_matrix_row (pnode->D, j);
00632 
00633     /* copy vector */
00634     gsl_vector_memcpy (v, &(u.vector));
00635 
00636     /* update parameters */
00637     gnn_node_local_update (node);
00638 
00639     return 0;
00640 }
00641 
00642 
00643 

Generated on Sun Jun 13 20:50:12 2004 for libgnn Gradient Retropropagation Machine Library by doxygen1.2.18