NFFT Logo 3.1.1 API Reference

nfsoft.c

00001 /*
00002  * Copyright (c) 2002, 2009 Jens Keiner, Stefan Kunis, Daniel Potts
00003  *
00004  * This program is free software; you can redistribute it and/or modify it under
00005  * the terms of the GNU General Public License as published by the Free Software
00006  * Foundation; either version 2 of the License, or (at your option) any later
00007  * version.
00008  *
00009  * This program is distributed in the hope that it will be useful, but WITHOUT
00010  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00011  * FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
00012  * details.
00013  *
00014  * You should have received a copy of the GNU General Public License along with
00015  * this program; if not, write to the Free Software Foundation, Inc., 51
00016  * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
00017  */
00018 
00019 /* $Id: nfsoft.c 3198 2009-05-27 14:16:50Z keiner $ */
00020 
00021 #include <stdio.h>
00022 #include <math.h>
00023 #include <string.h>
00024 #include <stdlib.h>
00025 #include <complex.h>
00026 #include "nfft3.h"
00027 #include "nfft3util.h"
00028 #include "infft.h"
00029 #include "wigner.h"
00030 
00031 #define DEFAULT_NFFT_CUTOFF    6
00032 #define FPT_THRESHOLD          1000
00033 
00034 void nfsoft_init(nfsoft_plan *plan, int N, int M)
00035 {
00036   nfsoft_init_advanced(plan, N, M, NFSOFT_MALLOC_X | NFSOFT_MALLOC_F
00037       | NFSOFT_MALLOC_F_HAT);
00038 }
00039 
00040 void nfsoft_init_advanced(nfsoft_plan *plan, int N, int M,
00041     unsigned int nfsoft_flags)
00042 {
00043   nfsoft_init_guru(plan, N, M, nfsoft_flags, PRE_PHI_HUT | PRE_PSI | MALLOC_X
00044       | MALLOC_F_HAT | MALLOC_F | FFTW_INIT | FFT_OUT_OF_PLACE,
00045       DEFAULT_NFFT_CUTOFF, FPT_THRESHOLD);
00046 }
00047 
00048 void nfsoft_init_guru(nfsoft_plan *plan, int B, int M,
00049     unsigned int nfsoft_flags, unsigned int nfft_flags, int nfft_cutoff,
00050     int fpt_kappa)
00051 {
00052   int N[3];
00053   int n[3];
00054 
00055   N[0] = 2* B + 2;
00056   N[1] = 2* B + 2;
00057   N[2] = 2* B + 2;
00058 
00059   n[0] = 8* B ;
00060   n[1] = 8* B ;
00061   n[2] = 8* B ;
00062 
00063   nfft_init_guru(&plan->nfft_plan, 3, N, M, n, nfft_cutoff, nfft_flags,
00064       FFTW_ESTIMATE | FFTW_DESTROY_INPUT);
00065 
00066   if ((plan->nfft_plan).nfft_flags & PRE_LIN_PSI)
00067   {
00068     nfft_precompute_lin_psi(&(plan->nfft_plan));
00069   }
00070 
00071   plan->N_total = B;
00072   plan->M_total = M;
00073   plan->fpt_kappa = fpt_kappa;
00074   plan->flags = nfsoft_flags;
00075 
00076   if (plan->flags & NFSOFT_MALLOC_F_HAT)
00077   {
00078     plan->f_hat = (C*) nfft_malloc((B + 1) * (4* (B +1)*(B+1)-1)/3*sizeof(C));
00079     if (plan->f_hat == NULL ) printf("Allocation failed!\n");
00080   }
00081 
00082   if (plan->flags & NFSOFT_MALLOC_X)
00083   {
00084     plan->x = (R*) nfft_malloc(plan->M_total*3*sizeof(R));
00085     if (plan->x == NULL ) printf("Allocation failed!\n");
00086   }
00087   if (plan->flags & NFSOFT_MALLOC_F)
00088   {
00089     plan->f = (C*) nfft_malloc(plan->M_total*sizeof(C));
00090       if (plan->f == NULL ) printf("Allocation failed!\n");
00091   }
00092 
00093   plan->wig_coeffs = (C*) nfft_malloc((nfft_next_power_of_2(B)+1)*sizeof(C));
00094   plan->cheby = (C*) nfft_malloc((2*B+2)*sizeof(C));
00095   plan->aux = (C*) nfft_malloc((2*B+4)*sizeof(C));
00096 
00097   if (plan->wig_coeffs == NULL ) printf("Allocation failed!\n");
00098   if (plan->cheby == NULL ) printf("Allocation failed!\n");
00099   if (plan->aux == NULL ) printf("Allocation failed!\n");
00100 
00101   plan->mv_trafo = (void (*) (void* ))nfsoft_trafo;
00102   plan->mv_adjoint = (void (*) (void* ))nfsoft_adjoint;
00103 
00104   plan->fpt_set = 0;
00105 }
00106 
00107 static void c2e(nfsoft_plan *my_plan, int even)
00108 {
00109   int j, N;
00110 
00112   N = 2* (my_plan ->N_total+1);
00113 
00115   my_plan->cheby[my_plan->N_total+1] = my_plan->wig_coeffs[0];
00116   my_plan->cheby[0]=0.0;
00117 
00118   for (j=1;j<my_plan->N_total+1;j++)
00119   {
00120     my_plan->cheby[my_plan->N_total+1+j]=0.5* my_plan->wig_coeffs[j];
00121     my_plan->cheby[my_plan->N_total+1-j]=0.5* my_plan->wig_coeffs[j];
00122   }
00123 
00124   C *aux= (C*) nfft_malloc((N+2)*sizeof(C));
00125 
00126   for(j=1;j<N;j++)
00127   aux[j]=my_plan->cheby[j];
00128 
00129   aux[0]=0.;
00130   aux[N]=0.;
00131 
00132   if (even>0)
00133   {
00134     my_plan->cheby[0]=(C) (-1.)/(2.*_Complex_I) * aux[1];
00135     for (j=1;j<N;j++)
00136     {
00137       my_plan->cheby[j]=(1./(2.*_Complex_I)*(aux[j+1]-aux[j-1]));
00138     }
00139 
00140   }
00141   free(aux);
00142   aux = NULL;
00143 }
00144 
00145 
00146 static fpt_set SO3_fpt_init(int l, fpt_set set, unsigned int flags, int kappa)
00147 {
00148   int N, t, k_start, k_end, k, m;
00149   int glo = 0;
00150   R *alpha, *beta, *gamma;
00151 
00153   if (flags & NFSOFT_USE_DPT)
00154   {
00155     if (l < 2)
00156       N = 2;
00157     else
00158       N = l;
00159 
00160     t = (int) log2(nfft_next_power_of_2(N));
00161 
00162   }
00163   else
00164   {
00166     if (l < 2)
00167       N = 2;
00168     else
00169       N = nfft_next_power_of_2(l);
00170 
00171     t = (int) log2(N);
00172   }
00173 
00175   alpha = (R*) nfft_malloc((N + 2) * sizeof(R));
00176   beta = (R*) nfft_malloc((N + 2) * sizeof(R));
00177   gamma = (R*) nfft_malloc((N + 2) * sizeof(R));
00178 
00180   if (flags & NFSOFT_NO_STABILIZATION)
00181   {
00182     set = fpt_init((2* N + 1) * (2* N + 1), t, 0U | FPT_NO_STABILIZATION);
00183   }
00184   else
00185   {
00186     set = fpt_init((2* N + 1) * (2* N + 1), t, 0U);
00187   }
00188 
00189   for (k = -N; k <= N; k++)
00190     for (m = -N; m <= N; m++)
00191     {
00193       k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
00194       k_end = N;
00195 
00196       SO3_alpha_row(alpha, N, k, m);
00197       SO3_beta_row(beta, N, k, m);
00198       SO3_gamma_row(gamma, N, k, m);
00199 
00200       fpt_precompute(set, glo, alpha, beta, gamma, k_start, kappa);
00201       glo++;
00202     }
00203 
00204   free(alpha);
00205   free(beta);
00206   free(gamma);
00207   alpha = NULL;
00208   beta = NULL;
00209   gamma = NULL;
00210 
00211   return set;
00212 }
00213 
00214 fpt_set SO3_single_fpt_init(int l, int k, int m, unsigned int flags, int kappa)
00215 {
00216   int N, t, k_start, k_end;
00217   R *alpha, *beta, *gamma;
00218   fpt_set set = 0;
00219 
00221   if (flags & NFSOFT_USE_DPT)
00222   {
00223     if (l < 2)
00224       N = 2;
00225     else
00226       N = l;
00227 
00228     t = (int) log2(nfft_next_power_of_2(N));
00229 
00230   }
00231   else
00232   {
00234     if (l < 2)
00235       N = 2;
00236     else
00237       N = nfft_next_power_of_2(l);
00238 
00239     t = (int) log2(N);
00240   }
00241 
00243   alpha = (R*) nfft_malloc((N + 2) * sizeof(R));
00244   beta = (R*) nfft_malloc((N + 2) * sizeof(R));
00245   gamma = (R*) nfft_malloc((N + 2) * sizeof(R));
00246 
00248   if (flags & NFSOFT_NO_STABILIZATION)
00249   {
00250     set = fpt_init(1, t, 0U | FPT_NO_STABILIZATION);
00251   }
00252   else
00253   {
00254     set = fpt_init(1, t, 0U);
00255   }
00256 
00258   k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
00259   k_end = N;
00260 
00261   SO3_alpha_row(alpha, N, k, m);
00262   SO3_beta_row(beta, N, k, m);
00263   SO3_gamma_row(gamma, N, k, m);
00264 
00265   fpt_precompute(set, 0, alpha, beta, gamma, k_start, kappa);
00266 
00267   free(alpha);
00268   free(beta);
00269   free(gamma);
00270   alpha = NULL;
00271   beta = NULL;
00272   gamma = NULL;
00273 
00274   return set;
00275 }
00276 
00277 void SO3_fpt(C *coeffs, fpt_set set, int l, int k, int m, unsigned int flags)
00278 {
00279   int N;
00281   C* x;
00283   C* y;
00284 
00285   int trafo_nr; 
00286   int k_start, k_end, j;
00287   int function_values = 0;
00288 
00290   if (flags & NFSOFT_USE_DPT)
00291   {
00292     N = l;
00293     if (l < 2)
00294       N = 2;
00295   }
00296   else
00297   {
00298     if (l < 2)
00299       N = 2;
00300     else
00301       N = nfft_next_power_of_2(l);
00302   }
00303 
00305   k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
00306   k_end = N;
00307   trafo_nr = (N + k) * (2* N + 1) + (m + N);
00308 
00310   x = (C*) nfft_malloc((k_end + 1) * sizeof(C));
00311 
00312   for (j = 0; j <= k_end; j++)
00313    x[j] = K(0.0);
00314 
00315 
00316   for (j = 0; j <= l - k_start; j++)
00317   {
00318     x[j + k_start] = coeffs[j];
00319   }
00320   for (j = l - k_start + 1; j <= k_end - k_start; j++)
00321   {
00322     x[j + k_start] = K(0.0);
00323   }
00324 
00326   y = (C*) nfft_malloc((k_end + 1) * sizeof(C));
00327 
00328   if (flags & NFSOFT_USE_DPT)
00329   { 
00330     dpt_trafo(set, trafo_nr, &x[k_start], y, k_end, 0U
00331         | (function_values ? FPT_FUNCTION_VALUES : 0U));
00332   }
00333   else
00334   { 
00335     fpt_trafo(set, trafo_nr, &x[k_start], y, k_end, 0U
00336         | (function_values ? FPT_FUNCTION_VALUES : 0U));
00337   }
00338 
00340   for (j = 0; j <= l; j++)
00341   {
00342     coeffs[j] = y[j];
00343   }
00344 
00347   free(x);
00348   free(y);
00349   x = NULL;
00350   y = NULL;
00351 }
00352 
00353 void SO3_fpt_transposed(C *coeffs, fpt_set set, int l, int k, int m,
00354     unsigned int flags)
00355 {
00356   int N, k_start, k_end, j;
00357   int trafo_nr; 
00358   int function_values = 0;
00360   C* x;
00362   C* y;
00363 
00366   if (flags & NFSOFT_USE_DPT)
00367   {
00368     N = l;
00369     if (l < 2)
00370       N = 2;
00371   }
00372   else
00373   {
00374     if (l < 2)
00375       N = 2;
00376     else
00377       N = nfft_next_power_of_2(l);
00378   }
00379 
00381   k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
00382   k_end = N;
00383   trafo_nr = (N + k) * (2* N + 1) + (m + N);
00384 
00386   y = (C*) nfft_malloc((k_end + 1) * sizeof(C));
00388   x = (C*) nfft_malloc((k_end + 1) * sizeof(C));
00389 
00390   for (j = 0; j <= l; j++)
00391   {
00392     y[j] = coeffs[j];
00393   }
00394   for (j = l + 1; j <= k_end; j++)
00395   {
00396     y[j] = K(0.0);
00397   }
00398 
00399   if (flags & NFSOFT_USE_DPT)
00400   {
00401     dpt_transposed(set, trafo_nr, &x[k_start], y, k_end, 0U
00402         | (function_values ? FPT_FUNCTION_VALUES : 0U));
00403   }
00404   else
00405   {
00406     fpt_transposed(set, trafo_nr, &x[k_start], y, k_end, 0U
00407         | (function_values ? FPT_FUNCTION_VALUES : 0U));
00408   }
00409 
00410   for (j = 0; j <= l; j++)
00411   {
00412     coeffs[j] = x[j];
00413   }
00414 
00416   free(x);
00417   free(y);
00418   x = NULL;
00419   y = NULL;
00420 }
00421 
00422 void nfsoft_precompute(nfsoft_plan *plan3D)
00423 {
00424   int j;
00425   int N = plan3D->N_total;
00426   int M = plan3D->M_total;
00427 
00430   for (j = 0; j < M; j++)
00431   {
00432     plan3D->nfft_plan.x[3* j ] = plan3D->x[3* j + 2];
00433     plan3D->nfft_plan.x[3* j + 1] = plan3D->x[3* j ];
00434     plan3D->nfft_plan.x[3* j + 2] = plan3D->x[3* j + 1];
00435   }
00436 
00437   for (j = 0; j < 3* plan3D ->nfft_plan.M_total; j++)
00438   {
00439     plan3D->nfft_plan.x[j] = plan3D->nfft_plan.x[j] * (1 / (2* PI ));
00440   }
00441 
00442   if ((plan3D->nfft_plan).nfft_flags & FG_PSI)
00443   {
00444     nfft_precompute_one_psi(&(plan3D->nfft_plan));
00445   }
00446   if ((plan3D->nfft_plan).nfft_flags & PRE_PSI)
00447   {
00448     nfft_precompute_one_psi(&(plan3D->nfft_plan));
00449   }
00450 
00452   plan3D->fpt_set = SO3_fpt_init(N, plan3D->fpt_set, plan3D->flags,
00453       plan3D->fpt_kappa);
00454 
00455   if ((plan3D->nfft_plan).nfft_flags & MALLOC_F_HAT)
00456   for (j = 0; j < plan3D->nfft_plan.N_total; j++)
00457     plan3D->nfft_plan.f_hat[j] = 0.0;
00458 
00459 }
00460 
00461 void nfsoft_trafo(nfsoft_plan *plan3D)
00462 {
00463   int i, j, m, k, max, glo1, glo2;
00464 
00465   i = 0;
00466   glo1 = 0;
00467   glo2 = 0;
00468 
00469   int N = plan3D->N_total;
00470   int M = plan3D->M_total;
00471 
00473   if (N == 0)
00474   {
00475     for (j = 0; j < M; j++)
00476       plan3D->f[j] = plan3D->f_hat[0];
00477     return;
00478   }
00479 
00480   for (j = 0; j < plan3D->nfft_plan.N_total; j++)
00481     plan3D->nfft_plan.f_hat[j] = 0.0;
00482 
00483   for (k = -N; k <= N; k++)
00484   {
00485     for (m = -N; m <= N; m++)
00486     {
00487 
00488       max = (ABS(m) > ABS(k) ? ABS(m) : ABS(k));
00489 
00490       for (j = 0; j <= N - max; j++)
00491       {
00492         plan3D->wig_coeffs[j] = plan3D->f_hat[glo1];
00493 
00494         if ((plan3D->flags & NFSOFT_NORMALIZED))
00495         {
00496           plan3D->wig_coeffs[j] = plan3D->wig_coeffs[j] * (1. / (2. * PI))
00497               * SQRT(0.5 * (2. * (max + j) + 1.));
00498         }
00499 
00500         if ((plan3D->flags & NFSOFT_REPRESENT))
00501         {
00502           if ((k < 0) && (k % 2))
00503           {
00504             plan3D->wig_coeffs[j] = plan3D->wig_coeffs[j] * (-1);
00505           }
00506           if ((m < 0) && (m % 2))
00507             plan3D->wig_coeffs[j] = plan3D->wig_coeffs[j] * (-1);
00508         }
00509 
00510         glo1++;
00511       }
00512 
00513       for (j = N - max + 1; j < nfft_next_power_of_2(N) + 1; j++)
00514         plan3D->wig_coeffs[j] = 0.0;
00515       //fprintf(stdout,"\n k= %d, m= %d \n",k,m);
00516       SO3_fpt(plan3D->wig_coeffs, plan3D->fpt_set, N, k, m, plan3D->flags);
00517 
00518       c2e(plan3D, ABS((k + m) % 2));
00519 
00520       for (i = 1; i <= 2* plan3D ->N_total + 2; i++)
00521       {
00522         plan3D->nfft_plan.f_hat[NFSOFT_INDEX(k, m, i - N - 1, N) - 1]
00523             = plan3D->cheby[i - 1];
00524         //fprintf(stdout,"%f \t", plan3D->nfft_plan.f_hat[NFSOFT_INDEX(k,m,i-N-1,N)-1]);
00525         //fprintf(stdout,"another index: %d for k=%d,m=%d,l=%d,N=%d \n", NFSOFT_INDEX(k,m,i-N-1,N)-1,k,m,i-N-1,N);
00526       }
00527 
00528     }
00529   }
00530 
00531   if (plan3D->flags & NFSOFT_USE_NDFT)
00532   {
00533     ndft_trafo(&(plan3D->nfft_plan));
00534   }
00535   else
00536   {
00537     nfft_trafo(&(plan3D->nfft_plan));
00538   }
00539 
00540   for (j = 0; j < plan3D->M_total; j++)
00541     plan3D->f[j] = plan3D->nfft_plan.f[j];
00542 
00543 }
00544 
00545 static void e2c(nfsoft_plan *my_plan, int even)
00546 {
00547   int N;
00548   int j;
00549 
00551   N = 2* (my_plan ->N_total+1);
00552   //nfft_vpr_complex(my_plan->cheby,N+1,"chebychev");
00553 
00554 
00555       if (even>0)
00556       {
00557         //my_plan->aux[N-1]= -1/(2*I)* my_plan->cheby[N-2];
00558         my_plan->aux[0]= 1/(2*_Complex_I)*my_plan->cheby[1];
00559 
00560         for(j=1;j<N-1;j++)
00561         {
00562           my_plan->aux[j]=1/(2*_Complex_I)*(my_plan->cheby[j+1]-my_plan->cheby[j-1]);
00563 }
00564 my_plan->aux[N-1]=1/(2*_Complex_I)*(-my_plan->cheby[j-1]);
00565 
00566 
00567 for(j=0;j<N;j++)
00568 {
00569 my_plan->cheby[j]= my_plan->aux[j];
00570 }
00571 }
00572 
00573 my_plan->wig_coeffs[0]=my_plan->cheby[my_plan->N_total+1];
00574 
00575 for(j=1;j<=my_plan->N_total;j++)
00576 {
00577 my_plan->wig_coeffs[j]=0.5*(my_plan->cheby[my_plan->N_total+j+1]+my_plan->cheby[my_plan->N_total+1-j]);
00578 }
00579 
00580 
00581 
00582 //nfft_vpr_complex(my_plan->wig_coeffs,my_plan->N_total,"chebychev ");
00583 
00584 }
00585 
00586 void nfsoft_adjoint(nfsoft_plan *plan3D)
00587 {
00588   int i, j, m, k, max, glo1, glo2;
00589 
00590   i = 0;
00591   glo1 = 0;
00592   glo2 = 0;
00593 
00594   int N = plan3D->N_total;
00595   int M = plan3D->M_total;
00596 
00597   //nothing much to be done for polynomial degree 0
00598   if (N == 0)
00599   {
00600     plan3D->f_hat[0]=0;
00601     for (j = 0; j < M; j++)
00602       plan3D->f_hat[0] += plan3D->f[j];
00603     return;
00604   }
00605 
00606   for (j = 0; j < M; j++)
00607   {
00608     plan3D->nfft_plan.f[j] = plan3D->f[j];
00609   }
00610 
00611   if (plan3D->flags & NFSOFT_USE_NDFT)
00612   {
00613     ndft_adjoint(&(plan3D->nfft_plan));
00614   }
00615   else
00616   {
00617     nfft_adjoint(&(plan3D->nfft_plan));
00618   }
00619 
00620   //nfft_vpr_complex(plan3D->nfft_plan.f_hat,plan3D->nfft_plan.N_total,"all results");
00621 
00622   glo1 = 0;
00623 
00624   for (k = -N; k <= N; k++)
00625   {
00626     for (m = -N; m <= N; m++)
00627     {
00628 
00629       max = (ABS(m) > ABS(k) ? ABS(m) : ABS(k));
00630 
00631       for (i = 1; i < 2* plan3D ->N_total + 3; i++)
00632       {
00633         plan3D->cheby[i - 1] = plan3D->nfft_plan.f_hat[NFSOFT_INDEX(k, m, i - N
00634             - 1, N) - 1];
00635       }
00636 
00637       //fprintf(stdout,"k=%d,m=%d \n",k,m);
00638       //nfft_vpr_complex(plan3D->cheby,2*plan3D->N_total+2,"euler");
00639       e2c(plan3D, ABS((k + m) % 2));
00640 
00641       //nfft_vpr_complex(plan3D->wig_coeffs,plan3D->N_total+1,"chebys");
00642       SO3_fpt_transposed(plan3D->wig_coeffs, plan3D->fpt_set, N, k, m,
00643           plan3D->flags);
00644       //nfft_vpr_complex(plan3D->wig_coeffs,plan3D->N_total+1,"wigners");
00645       //  SO3_fpt_transposed(plan3D->wig_coeffs,N,k,m,plan3D->flags,plan3D->fpt_kappa);
00646 
00647 
00648       for (j = max; j <= N; j++)
00649       {
00650         if ((plan3D->flags & NFSOFT_REPRESENT))
00651         {
00652           if ((k < 0) && (k % 2))
00653           {
00654             plan3D->wig_coeffs[j] = -plan3D->wig_coeffs[j];
00655           }
00656           if ((m < 0) && (m % 2))
00657             plan3D->wig_coeffs[j] = -plan3D->wig_coeffs[j];
00658         }
00659 
00660         plan3D->f_hat[glo1] = plan3D->wig_coeffs[j];
00661 
00662         if ((plan3D->flags & NFSOFT_NORMALIZED))
00663         {
00664           plan3D->f_hat[glo1] = plan3D->f_hat[glo1] * (1 / (2. * PI)) * SQRT(
00665               0.5 * (2. * (j) + 1.));
00666         }
00667 
00668         glo1++;
00669       }
00670 
00671     }
00672   }
00673 }
00674 
00675 void nfsoft_finalize(nfsoft_plan *plan)
00676 {
00677   /* Finalise the nfft plan. */
00678   nfft_finalize(&plan->nfft_plan);
00679   free(plan->wig_coeffs);
00680   free(plan->cheby);
00681   free(plan->aux);
00682 
00683   fpt_finalize(plan->fpt_set);
00684   plan->fpt_set = NULL;
00685 
00686   if (plan->flags & NFSOFT_MALLOC_F_HAT)
00687   {
00688     //fprintf(stderr,"deallocating f_hat\n");
00689     free(plan->f_hat);
00690   }
00691 
00692   /* De-allocate memory for samples, if neccesary. */
00693   if (plan->flags & NFSOFT_MALLOC_F)
00694   {
00695     //fprintf(stderr,"deallocating f\n");
00696     free(plan->f);
00697   }
00698 
00699   /* De-allocate memory for nodes, if neccesary. */
00700   if (plan->flags & NFSOFT_MALLOC_X)
00701   {
00702     //fprintf(stderr,"deallocating x\n");
00703     free(plan->x);
00704   }
00705 }
00706 
00707 int posN(int n, int m, int B)
00708 {
00709   int pos;
00710 
00711   if (n > -B)
00712     pos = posN(n - 1, m, B) + B + 1 - MAX(ABS(m), ABS(n - 1));
00713   else
00714     pos = 0;
00715   //(n > -B? pos=posN(n-1,m,B)+B+1-MAX(ABS(m),ABS(n-1)): pos= 0)
00716   return pos;
00717 }
00718 

Generated on 17 Aug 2009 by Doxygen 1.5.3