source: GLBFrejus/HEAD/src/glb_minimize.c @ 157

Last change on this file since 157 was 2, checked in by campagne, 19 years ago

first import

File size: 28.7 KB
Line 
1/* GLoBES -- General LOng Baseline Experiment Simulator
2 * (C) 2002 - 2004,  The GLoBES Team
3 *
4 * GLoBES is mainly intended for academic purposes. Proper
5 * credit must be given if you use GLoBES or parts of it. Please
6 * read the section 'Credit' in the README file.
7 *
8 * This program is free software; you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License as published by
10 * the Free Software Foundation; either version 2 of the License, or
11 * (at your option) any later version.
12 *
13 * This program is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public License
19 * along with this program; if not, write to the Free Software
20 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
21 */
22
23
24
25
26
27
28/* --------------------------------------------------------
29 * --- Comment on numerical accuracy ----------------------
30 * --------------------------------------------------------
31 * The whole code uses double precision, whatever this is
32 * on Your machine. The following information was derived
33 * on 32-bit Pentium III processor.
34 * Even though double precision should give relative errors of
35 * order <1E-12, functions involving minimization like ChiTheta
36 * may be off in the order of 1E-8 (absolute). This may cause
37 * ChiTheta and SingleChiTheta to yield results differing by 1E-8.
38 * The same applies for all other ChiglbXSection and SingleChiglbXSection functions.
39 * The crucial lines are:
40 * erg2 = erg2
41 *    + glb_prior(x[1],start[2],inp_errs[2])
42 *    + glb_prior(x[3],start[4],inp_errs[4])
43 *    + glb_prior(x[4],start[5],inp_errs[5])
44 *    + glb_prior(x[5],(glb_experiment_list[glb_single_experiment_number]).density_center,
45 *      (glb_experiment_list[glb_single_experiment_number]).density_error);
46 * in the chi_xxx functions. This looks of course different
47 * in the MDchi_xxx functions. Basically its a matter of how
48 * and in which order the truncations are performed. This may
49 * also change under compiler optimization. In any case the errors
50 * should be very small (i.e. <<1E-6) and should not have any impact
51 * on the physics results.
52 */
53
54
55#include <stdio.h>
56#include <stdlib.h>
57#include <math.h>
58#include <setjmp.h>
59#include <globes/globes.h>
60
61#include "glb_probability.h"
62#include "glb_fluxes.h"
63#include "glb_rate_engine.h"
64#include "glb_min_sup.h"
65#include "glb_types.h"
66#include "glb_multiex.h"
67#include "glb_error.h"
68#include "glb_wrapper.h"
69
70
71
72
73
74/* Warning-- fiddling around with these in general may
75 * influence the stability of the minimization process !
76 *
77 * The latest Recator (D-CHOOZ as in hep-ph/0403068) setups
78 * have TOLSYS 1.0E-12, whereas all the other setups have TOLSYS 1.0E-5
79 * and TOLSYS 1.0E-8 may be faster by 20% for Nufact (improves convergence
80 * of the second layer of minimization)
81 */
82
83#define TOLSYS 1.0E-8
84#define TOLOSC 1.0E-5
85
86
87/* defines wether the program will run
88 * with (1) or w/o mathematica (0)
89 */
90
91#define GLB_MATHLINK 1
92
93/* Defines wether the thing runs on UNIX (0) or not (1) */
94
95#define UNIX 0
96
97#define FLOAT double
98
99//----------------------
100
101/* this is for making the Chi... functions abortable */
102static int okay_flag = 0;
103/* this needed by the wrappers in order not to return garbage */
104static jmp_buf env; 
105/* this is need by setjmp() longjmp() */
106//--------------------------------------------------------
107
108
109/* global variabels */
110int glb_single_experiment_number=0;
111
112
113static double inp_errs[7];
114static double start[7];
115static int count=0;
116static int errordim[32];
117static double sysglb_calc_buffer[6];
118
119//-------------------------------------------------------------
120//------------------    MLAbort handling     ------------------
121//-------------------------------------------------------------
122// the whole stuff is being switched of if GLB_MATHLINK is set to 0.
123
124// this function handels a caught abort message
125// and jumps to the point where setjmp(env) has
126// been called. usually from one of the ChiTheta(),
127// MDChiTheta() or Chi() or SingleChi(). this
128// function then will return a garbage value. this
129// fact is comunicated to the wrappers with the
130// okay_flag. if they encounter an okay_flag=1 they
131// will send the Abort[] mathematica comand and then
132// return.
133
134static void ml_abort_handler()
135{
136  fprintf(stderr,"Mathlink programe catched abort!\n");
137  longjmp(env,1);
138}
139
140// this function checks wether a abort message is sent by the
141// math-kernel, if so it calls the ml_abort_handler().
142// otherwise it calls MLCallYieldFunction() in order to give
143// the kernel a possibility to send its message (this should
144// be needed only in Windows and on MacOS)
145
146#ifdef MLabort
147static void ml_abort_check(int flag)
148{
149
150  if(flag==0) return;
151  if(!MLAbort)
152    {
153      if(UNIX==1) MLCallYieldFunction(MLYieldFunction(stdlink),
154                                      stdlink,(MLYieldParameters)0);
155    }
156  if(MLAbort)
157    {
158      ml_abort_handler();
159    }
160}
161#else
162static void  ml_abort_check(int flag)
163{
164  /* To silence gcc -Wall */
165  int i;
166  i=flag;
167  return;
168}
169
170#endif
171
172
173// functions for returning the bfp of systematic parameters
174
175static void to_glb_calc_buffer(double sys[],int len)
176{
177  int i;
178  for (i=0;i<len;i++) sysglb_calc_buffer[i]=sys[i+1];
179  return;
180}
181
182static void clear_glb_calc_buffer()
183{
184  int i;
185  for (i=0;i<6;i++) sysglb_calc_buffer[i]=0.0;
186  return;
187}
188
189static double* read_glb_calc_buffer()
190{
191  return &sysglb_calc_buffer[0];
192}
193
194//--------------------------------------------------------------------------
195//--------------------------------------------------------------------------
196//------------------------ Systematics -------------------------------------
197//--------------------------------------------------------------------------
198//--------------------------------------------------------------------------
199
200
201// this an init function -
202// setting up the input matrix for minimization
203// ie. the set of starting directions
204
205
206static void init_mat(double **m, int dim)
207{
208  int i;
209  int j;
210  for (i=1 ;i<= dim;i++)
211    {
212      for (j=1;j<= dim;j++)
213        {
214          if (i==j)
215            {
216              m[i][j]=1.0;
217            }
218          else
219            {
220              m[i][j]=0.0;
221            }   
222        }
223    }
224}
225
226//This is a wraper for accessing the minimizer
227
228static double chi_sys_wrap(double (*chi_func)(), int dimension, int alpha_in)
229{
230  double **mat; // contains the direction set
231  double *sp; // stores the coordinates of the minimum
232  double result=0; // stores the minimum value
233  int it=0; // counts the number of iterations
234  mat=glb_alloc_mat(1,dimension,1,dimension);
235  sp=glb_alloc_vec(1,6); // that is the maximal length
236  glb_rule_number=alpha_in;
237  init_mat(mat,dimension);
238
239  sp[1]=1;
240  sp[2]=0;
241  sp[3]=glb_bg_norm_center[alpha_in];
242  sp[4]=glb_bg_tilt_center[alpha_in];
243  sp[5]=glb_tre_null_center[alpha_in];
244  sp[6]=glb_tre_tilt_center[alpha_in];
245  glb_powell(sp,mat,dimension,TOLSYS,&it,&result,chi_func); 
246  glb_free_vec(sp,1,6);
247  glb_free_mat(mat,1,dimension,1,dimension);
248  return result;
249}
250//Chi^2 where the systematics are integrated out
251
252void glb_set_errordim(int typ, int rule)
253{
254  errordim[rule]=typ;
255}
256
257int glb_check_errordim(int rule)
258{
259  return errordim[rule];
260}
261
262
263
264// probably will be replaced by a more object oriented approach
265
266// those two are needed in order to write unique
267// chi^2 functions
268double *glb_sys_errors;
269double *glb_sys_centers;
270
271struct glb_systematic glb_init_systematic(double (*chi_func)(),int dimension,
272                            double *sp, double *errors, double (*evalf)(),
273                            char info[])
274{
275  struct glb_systematic out;
276  double *si,*se;
277  int i;
278  si=(double *) glb_malloc(sizeof(double)*dimension);
279  se=(double *) glb_malloc(sizeof(double)*dimension);
280  for(i=0;i<dimension;i++) si[i]=sp[i];
281  for(i=0;i<dimension;i++) se[i]=errors[i];
282  out.chi_func=chi_func;
283  out.dimension=dimension;
284  out.sp=si;
285  out.errors=se;
286  out.evalf=evalf;
287  out.info=info;
288  return out;
289}
290
291
292
293double glb_evaluate_chi(struct glb_systematic *in)
294{
295  int i;
296  double **mat; // contains the direction set
297  double *spi; // stores the coordinates of the minimum
298  double result=0; // stores the minimum value
299  int it=0; // counts the number of iterations
300  mat=glb_alloc_mat(1,in->dimension,1,in->dimension);
301  spi=glb_alloc_vec(1,in->dimension); // that is the maximal length
302  init_mat(mat,in->dimension);
303  for(i=1;i<in->dimension+1;i++) spi[i]=in->sp[i-1];
304  // this part of the unified interface to this quantities
305  glb_sys_centers=in->sp;
306  glb_sys_errors=in->errors;
307  //------------------------
308  glb_powell(spi,mat,in->dimension,TOLSYS,&it,&result,in->chi_func); 
309  glb_free_vec(spi,1,in->dimension);
310  glb_free_mat(mat,1,in->dimension,1,in->dimension);
311  return result;
312}
313
314// This is the central switch board for choosing the
315// chi^2 function and  systematics treatment according
316// to each value of errordim
317
318static double chi_dispatch(int error_dim, int i)
319{     
320  // finally
321  double erg2=0;
322 
323  if(error_dim==0) //chi^2 with 4 sys. parameters
324    {
325      erg2=chi_sys_wrap(&glb_chi_sys_w_bg,4,i);
326     
327    }
328  else if(error_dim==1) //chi^2 with 6 sys. parameters
329    {
330      erg2=chi_sys_wrap(&glb_chi_sys_w_bg2,6,i); 
331    }
332  else if(error_dim==2) // chi^2 without sys. parameters
333    {
334      erg2=chi_sys_wrap(&glb_chi_sys_w_bg,0,i); 
335      // it seems that the minimizer correctly interprest zero dimensions
336      // as the function value at the starting point 
337    }
338  else if(error_dim==3) // very special for JHF-HK
339    {
340      if(i<4) erg2=chi_sys_wrap(&glb_chi_sys_w_bg,4,i); 
341      else erg2=chi_sys_wrap(&glb_chi_sys_w_bgtot,4,i); 
342    }
343  else if(error_dim==4) // total chi^2 with 4 sys parameters
344    {
345      erg2=chi_sys_wrap(&glb_chi_sys_w_bgtot2,4,i); 
346    }
347  else if(error_dim==5) // obsolete
348    {
349      fprintf(stderr,"Warning: obsolete errordim\n");
350      erg2=0;
351    }
352  else if(error_dim==6) // obsolete
353    {
354      fprintf(stderr,"Warning: obsolete errordim\n");
355      erg2=0;
356    }
357  else if(error_dim==7) // free normalization, i.e. spectrum only
358    {
359      erg2=chi_sys_wrap(&glb_chi_spec,1,i); 
360    }
361  else if(error_dim==8) //total chi^2 w/o systematics
362    {
363      erg2=chi_sys_wrap(&glb_chi_sys_w_bgtot2,0,i); 
364    }
365  else if(error_dim==9) // chi^2 with 4 syst params and
366    // true energy calibration
367    {
368      erg2=chi_sys_wrap(&glb_chi_sys_w_bg_calib,4,i); 
369    }
370  else if(error_dim==10) // total chi^2 with 4 syst params
371    {
372      erg2=chi_sys_wrap(&glb_chi_sys_w_bgtot,4,i); 
373    }
374  else if(error_dim==20) // total chi^2 with 4 syst params
375    {
376      erg2=glb_evaluate_chi(&sys_calc[i]); 
377    }
378  else if(error_dim==21) // total chi^2 with 4 syst params
379    {
380      erg2=sys_calc[i].evalf(&sys_calc[i]); 
381    }
382  else
383    {
384      erg2=0;
385    }
386return erg2;
387}
388
389static double ChiS0(int typ[])
390{
391  double erg2;
392  int i,rul;
393 
394  ml_abort_check(GLB_MATHLINK);
395 
396  erg2=0;
397  for (i=0;i<glb_num_of_rules;i++) 
398    {
399      erg2+=chi_dispatch(typ[i],i);
400    }
401 
402  glb_rule_number=0;
403  return erg2;
404}
405
406// this is the same as ChiSO()  but it allows access to a single rule !
407
408static double ChiS0_Rule(int typ[],int rule)
409{
410  double erg2;
411  clear_glb_calc_buffer();
412  erg2=chi_dispatch(typ[rule],rule);
413  glb_rule_number=0;
414  return erg2;
415}
416
417//---------------------------------------------------------------
418//---------- MES begins here ------------------------------------
419//---------------------------------------------------------------
420
421// redfinition of ChiS
422
423static double ChiS()
424{
425  double erg;
426  int i;
427  erg=0;
428  for (i=0;i<glb_num_of_exps;i++)
429    {
430      glbSetExperiment(glb_experiment_list[i]);
431      erg += ChiS0(glb_experiment_list[i]->errordim);
432    }
433  return erg;
434}
435
436//redefinition of Chi
437
438static double Chi(double x[])
439{
440  int i; 
441  double erg;
442  glb_set_c_vacuum_parameters(x[0],x[1],x[2],x[3]);
443  glb_set_c_squared_masses(0,x[4],x[5]); 
444  for (i=0;i<glb_num_of_exps;i++)
445    {
446      glbSetExperiment(glb_experiment_list[i]);
447      glb_set_profile_scaling(x[6+i],i);
448      glb_set_new_rates();
449    }
450  if (setjmp(env)==1) 
451    {
452      okay_flag=1;
453      return erg; 
454    }
455  erg=ChiS();
456  return erg;
457}
458
459
460
461
462// chi^2 with systematics for each Experiment
463
464static double SingleChi(double x[7],int exp)
465{
466  double erg;
467  glb_set_c_vacuum_parameters(x[0], x[1],x[2],x[3]);
468  glb_set_c_squared_masses(0,x[4],x[5]);
469
470      glbSetExperiment(glb_experiment_list[exp]);
471      glb_set_profile_scaling(x[6],exp);
472      glb_set_new_rates();
473     
474
475
476  glbSetExperiment(glb_experiment_list[exp]);
477  if (setjmp(env)==1) 
478    {
479      okay_flag=1;
480      return erg; 
481    }
482  erg=ChiS0(errordim);
483  return erg;
484}
485
486// chi^2 with systematics for each Experiment
487
488static double SingleRuleChi(double x[7],int exp, int rule)
489{
490  double erg;
491  glb_set_c_vacuum_parameters(x[0], x[1],x[2],x[3]);
492  glb_set_c_squared_masses(0,x[4],x[5]);
493
494      glbSetExperiment(glb_experiment_list[exp]);
495      glb_set_profile_scaling(x[6],exp);
496      glb_set_new_rates();
497     
498
499
500  glbSetExperiment(glb_experiment_list[exp]);
501  erg=ChiS0_Rule(errordim,rule);
502  return erg;
503}
504
505/* Wrappers for the API */
506
507double glbChiSys(const glb_params in,int experiment, int rule)
508{
509  int i;
510  double res,x[38];
511 
512  if(in==NULL) 
513    {
514      glb_error("Failure in glbChiSys: Input pointer must be non-NULL");
515      return -1;
516    }
517  for(i=0;i<GLB_OSCP;i++) x[i]=glbGetOscParams(in,i);
518  if(experiment==GLB_ALL)
519    {
520      if(rule==GLB_ALL)
521        {
522          for(i=0;i<glb_num_of_exps;i++) 
523            x[i+GLB_OSCP]=glbGetDensityParams(in,i);
524          res=Chi(x);
525        }
526      else
527        {
528          res=0;
529           for(i=0;i<glb_num_of_exps;i++)
530             {   
531               if(rule < glb_experiment_list[i]->numofrules)
532                 {
533                   x[GLB_OSCP]=glbGetDensityParams(in,i);
534                   res+=SingleRuleChi(x,i,rule);
535                 }
536               else
537                 {glb_error("Invalid rule number");return -1;}
538             }
539        }
540    }
541  else
542    {
543      if(experiment >= glb_num_of_exps) 
544        {
545          glb_error("Failure in glbChiSys: 2nd argument must be smaller than"
546                    "glb_num_of_exps");
547          return -1;
548        }
549      x[GLB_OSCP]=glbGetDensityParams(in,experiment);
550      if(rule==GLB_ALL) res=SingleChi(x,experiment);
551        else
552          {
553            if(rule >= glb_experiment_list[experiment]->numofrules)
554              {
555                glb_error("Failure in glbChiSys: 3rd argument must be"
556                          " smaller than numofrules");
557                return -1;
558              }
559            res=SingleRuleChi(x,experiment,rule);
560          }
561    }
562  return res;
563}
564
565
566//-----------------------------------------------------------------
567//------------------- END -----------------------------------------
568//--------------- Systematics -------------------------------------
569//-----------------------------------------------------------------
570
571
572//--------------------------------------------------------------
573//--------- Setting the glb_priors ---------------------------------
574//--------------------------------------------------------------
575
576int glb_set_input_errors(double a, double b, double c, double d, double e, double f)
577{
578 
579  inp_errs[1]=a;
580  inp_errs[2]=b;
581  inp_errs[3]=c;
582  inp_errs[4]=d;
583  inp_errs[5]=e;
584  inp_errs[6]=f;
585  return 0;
586}
587
588int glb_set_starting_values(double a, double b, double c, double d, double e, double f)
589{
590  start[1]=a;
591  start[2]=b;
592  start[3]=c;
593  start[4]=d;
594  start[5]=e;
595  start[6]=f;
596  return 0;
597}
598
599
600
601
602
603
604// setting starting values for densities and errors on densitites
605// for different experiments
606
607void glbSetDensityPrior(double start, double error, int typ)
608{
609  glb_experiment_list[typ]->density_center=start;
610  glb_experiment_list[typ]->density_error=error;
611}
612
613void glbSetDensityStartingValue(double start, int typ)
614{
615  glb_experiment_list[typ]->density_center=start;
616}
617
618
619void glbSetDensityInputError(double error, int typ)
620{
621  glb_experiment_list[typ]->density_error=error;
622}
623
624
625double* glb_return_input_errors()
626{
627  double* out;
628  int i;
629  i=6+glb_num_of_exps;
630  out=(double*) glb_malloc(i*sizeof(double));
631  for(i=0;i<6;i++) out[i]=inp_errs[i];
632  for(i=0;i<glb_num_of_exps;i++) out[i+6]=glb_experiment_list[i]->density_error;
633  return out;
634}
635
636double* glb_return_input_values()
637{
638  double* out;
639  int i;
640  i=6+glb_num_of_exps;
641  out=(double*) glb_malloc(i*sizeof(double));
642  for(i=0;i<6;i++) out[i]=start[i];
643  for(i=0;i<glb_num_of_exps;i++) out[i+6]=glb_experiment_list[i]->density_center;
644  return out;
645}
646 
647 
648
649//--------------------------------------------------------------
650//--------- Setting the glb_priors for th12 ------------------------
651//--------------------------------------------------------------
652
653int glb_set_solar_input_errors(double a)
654{
655  inp_errs[0]=a;
656  return 0;
657}
658
659int glb_set_solar_starting_values(double a)
660{
661  start[0]=a;
662  return 0;
663}
664
665
666
667//----------------------------------------------------------------
668//----------- shifted rate access --------------------------------
669//----------------------------------------------------------------
670
671static void ReturnShiftedRates(int rulenumber, int exp, double* inrates)
672{
673  int i;
674  int bins;
675  glbSetExperiment(glb_experiment_list[exp]);
676  //glb_set_new_rates();
677  //glbSetExperiment(&glb_experiment_list[exp]);
678  //ChiS0_Rule(errordim,rulenumber);
679  bins=rulenumber;
680  bins=glb_get_number_of_bins();
681  for(i=0;i<bins;i++) inrates[i]=glb_chirate[i];
682
683}
684
685//------------------------------------------------------------------
686//---- Chi^2 with arbitrary number of free parameters --------------
687//----------------------- 23.01.2004 -------------------------------
688//------------------------------------------------------------------
689
690//This serves for ChiNP
691
692static double fix_params[37];
693static int para_tab[37];
694static int index_tab[37];
695static int n_free;
696static int n_fix;
697
698//-----------------------
699
700
701static void SelectProjection(int *vec)
702{
703  int i,c,c1,c2;
704
705  for(i=0;i<6+glb_num_of_exps;i++) para_tab[i]=vec[i];
706  c=0;
707  c1=0;
708  c2=0;
709  for(i=0;i<6+glb_num_of_exps;i++)
710    {
711      if(para_tab[i]==1) c++;
712    }
713  for(i=0;i<6+glb_num_of_exps;i++)
714    {
715      if(para_tab[i]==1) 
716        {
717          index_tab[c1]=i;
718          c1++;
719        }
720      else if(para_tab[i]==0)
721        {
722          index_tab[c+c2]=i;
723          c2++;
724        }
725      else
726        {
727          glb_fatal("SelectProjection input error\n");
728        } 
729    }
730  n_free=c;
731  n_fix=c2;
732  return;
733}
734
735
736static int CheckFree()
737{
738  int k;
739  k=n_free;
740  return k;
741}
742
743static int* CheckProjection()
744{
745 
746  return &para_tab[0];
747}
748
749
750static double sglb_prior(double x, double center, double sigma)
751{
752  if(fabs(sigma-0)<1E-12) return 0;
753  return (x-center)*(x-center)/sigma/sigma;
754} 
755
756// multi-experiment functions MDglbXSection
757static double MD_chi_NP(double x[])
758{
759  double erg2;
760  double y[37];
761  int i; 
762  count = count +1;
763  for(i=0;i<n_free;i++) y[index_tab[i]]=x[i+1]; 
764  // This basically is superflous, however it appears to be safer not
765  // change a global (i.e to this file) parameter (fix_params) at this place
766  for(i=n_free;i<n_free+n_fix;i++) y[index_tab[i]]=fix_params[index_tab[i]]; 
767
768  //fprintf(stderr,"x2 %f\n",x2[1]);
769  glb_set_c_vacuum_parameters(y[0],y[1],y[2],y[3]);
770  glb_set_c_squared_masses(0,y[4],y[5]);
771  for (i=0;i<glb_num_of_exps;i++)
772    {
773      glbSetExperiment(glb_experiment_list[i]);
774      glb_set_profile_scaling(y[6+i],i);
775      glb_set_new_rates();
776    }
777 
778  erg2=ChiS();
779 
780  // adding the glb_priors
781  for(i=0;i<n_free;i++)
782    {
783      if(index_tab[i]<6)
784        {
785          erg2+=sglb_prior(y[index_tab[i]],
786                      start[index_tab[i]],inp_errs[index_tab[i]]); 
787         
788        }
789      else
790        {
791          erg2+=sglb_prior(y[index_tab[i]],
792                     (glb_experiment_list[index_tab[i]-6])->density_center,
793               (glb_experiment_list[index_tab[i]-6])->density_error);
794        }
795    }
796 
797  return erg2;
798}
799 
800
801// single-experiment functions ChiglbXSection
802//This serves for ChiNP
803
804static double s_fix_params[7];
805static int s_para_tab[7];
806static int s_index_tab[7];
807static int s_n_free;
808static int s_n_fix;
809
810//-----------------------
811
812
813static void single_SelectProjection(int set)
814{
815  int i,c,c1,c2;
816
817  for(i=0;i<6;i++) s_para_tab[i]=para_tab[i];
818  s_para_tab[6]=para_tab[6+set];
819  c=0;
820  c1=0;
821  c2=0;
822  for(i=0;i<6+1;i++)
823    {
824      if(s_para_tab[i]==1) c++;
825    }
826  for(i=0;i<6+1;i++)
827    {
828      if(s_para_tab[i]==1) 
829        {
830          s_index_tab[c1]=i;
831          c1++;
832        }
833      else if(s_para_tab[i]==0)
834        {
835          s_index_tab[c+c2]=i;
836          c2++;
837        }
838      else
839        {
840          fprintf(stderr,"SelectProjection input error\n");
841        } 
842    }
843  s_n_free=c;
844  s_n_fix=c2;
845  return;
846}
847
848//------------------ JEC START: 13/5/05 ----------------
849int same_th23(double start_th23, double fit_th23) {
850 
851  // testing if fit_th23 is at the same side as the true value
852 
853  return  ( 
854           (start_th23 > M_PI/4. && fit_th23 > M_PI/4.) ||
855           (start_th23 < M_PI/4. && fit_th23 < M_PI/4.) 
856           );
857}//same_th23
858
859//---------------------------------------------------------------------
860int same_hier(double start_dmq, double fit_dmq) {
861 
862  // testing if fit_dmq has the same sign as the true value
863  return ( 
864          (fit_dmq * start_dmq) > 0 
865          );
866}//same_hier
867
868//---------------------------------------------------------------------
869
870/* bool test(double start_dmq, double fit_dmq,  */
871/*        double start_th23, double fit_th23) {     */
872/*   //test if Hierarchy is the same, and if Octant is the same */
873/*   return same_th23(start_th23,fit_th23) && same_hier(start_dmq,fit_dmq); */
874/* } */
875//------------------ JEC END: 13/5/05 ----------------
876
877
878static double chi_NP(double x[])
879{
880  double erg2;
881  double y[7];
882  int i; 
883  count = count +1;
884  for(i=0;i<s_n_free;i++) y[s_index_tab[i]]=x[i+1]; 
885  // This basically is superflous, however it appears to be safer not
886  // change a global (i.e to this file) parameter (fix_params) at this place
887  for(i=s_n_free;i<s_n_free+s_n_fix;i++) y[s_index_tab[i]]
888                                         =s_fix_params[s_index_tab[i]]; 
889
890  //fprintf(stderr,"x2 %f\n",x2[1]);
891  glb_set_c_vacuum_parameters(y[0],y[1],y[2],y[3]);
892  glb_set_c_squared_masses(0,y[4],y[5]);
893 
894  glbSetExperiment(glb_experiment_list[glb_single_experiment_number]);
895  glb_set_profile_scaling(y[6],glb_single_experiment_number);
896  glb_set_new_rates();
897   
898 
899  erg2=ChiS0(errordim);
900 
901 
902  // adding the glb_priors
903  for(i=0;i<s_n_free;i++)
904    {
905      if(s_index_tab[i]<6)
906        {
907          erg2+=sglb_prior(y[s_index_tab[i]],
908                      start[s_index_tab[i]],inp_errs[s_index_tab[i]]); 
909          //----- JEC START: 13/5/05 ----
910          if(2 == s_index_tab[i]) { //Theta23 special case
911            if ( !same_th23(start[s_index_tab[i]],y[s_index_tab[i]]) ) {
912              erg2 += 1.e6;
913            }
914          } else if (5 == s_index_tab[i]) {
915            if ( !same_hier(start[s_index_tab[i]],y[s_index_tab[i]]) ) {
916              erg2 += 1.e6;
917            }
918          }
919          //----- JEC END: 13/5/05 ----
920        }
921      else
922        {
923          erg2+=sglb_prior(y[6],
924                     (glb_experiment_list[glb_single_experiment_number])->density_center,
925               (glb_experiment_list[glb_single_experiment_number])->density_error);
926        }
927    }
928 
929  return erg2;
930}
931 
932
933
934/* This implemenst the API for the ChiXXX functions */
935
936
937 
938static double internal_glbSingleChiNP(const glb_params in, glb_params out, 
939                               int exp)
940{
941  double *sp2;
942  double **mat2;
943  double er1;
944  glb_projection fbuf,fnew; 
945
946  double x[38];
947  int it;
948  int i;
949  int dim;
950  fbuf=glbAllocProjection();
951  fnew=glbAllocProjection();
952  if(exp >=  glb_num_of_exps)
953    {
954      glb_error("Failure in internal_glbSingleChiNP: exp must be smaller than"
955                " glb_num_of_exps");
956      return -1;
957    }
958
959 
960
961  glb_single_experiment_number=exp;
962 
963  //creating memory
964  single_SelectProjection(exp);
965
966  dim=s_n_free;
967  /* declaring temporariliy all densities of all other experiments as fixed */
968  glbGetProjection(fbuf);
969  fnew=glbCopyProjection(fbuf,fnew);
970  fnew=glbSetDensityProjectionFlag(fnew,GLB_FIXED,GLB_ALL);
971  fnew=glbSetDensityProjectionFlag(fnew,glbGetDensityProjectionFlag(fbuf,exp)
972                              ,exp); 
973  glbSetProjection(fnew); 
974  /* - finis - */
975  mat2=glb_alloc_mat(1,dim,1,dim);
976  sp2=glb_alloc_vec(1,dim);
977  init_mat(mat2,dim);
978  //initializing various things
979  count=0;
980 
981  if(out!=NULL) {
982    out=glbCopyParams(in,out);
983    if(out==NULL) 
984      {
985        glb_error("Failure while copying input of glbChiNP");
986        return -1;
987      }
988  }
989 
990  for(i=0;i<GLB_OSCP;i++) x[i]=glbGetOscParams(in,i);
991  x[GLB_OSCP]=glbGetDensityParams(in,exp);
992
993  for(i=0;i<GLB_OSCP+1;i++) s_fix_params[i]=x[i];
994  for(i=0;i<s_n_free;i++) sp2[i+1]=x[s_index_tab[i]];
995  if (setjmp(env)==1) 
996    {
997      okay_flag=1;
998      return er1; 
999    }
1000  glb_powell2(sp2,mat2,dim,TOLOSC,&it,&er1,chi_NP);
1001  if(out!=NULL)
1002    {
1003      for(i=0;i<s_n_free;i++) 
1004        {
1005          if(s_index_tab[i]<GLB_OSCP)
1006            glbSetOscParams(out,sp2[i+1],s_index_tab[i]);
1007          else
1008            glbSetDensityParams(out,sp2[i+1],exp);
1009        }
1010      out=glbSetIteration(out,count);
1011    }
1012 
1013  glb_free_vec(sp2,1,dim);
1014  glb_free_mat(mat2,1,dim,1,dim);
1015  glbSetProjection(fbuf);
1016  glbFreeProjection(fnew);
1017  glbFreeProjection(fbuf);
1018  return er1;
1019}
1020 
1021
1022static double internal_glbChiNP(const glb_params in, glb_params out)
1023{
1024  double *sp2;
1025  double **mat2;
1026  double er1;
1027 
1028
1029  double x[38];
1030  int it;
1031  int i;
1032  int dim;
1033  //creating memory
1034 
1035 
1036  dim=n_free;
1037 
1038  mat2=glb_alloc_mat(1,dim,1,dim);
1039  sp2=glb_alloc_vec(1,dim);
1040  init_mat(mat2,dim);
1041  //initializing various things
1042  count=0;
1043 
1044
1045 if(out!=NULL) {
1046    out=glbCopyParams(in,out);
1047    if(out==NULL) 
1048      {
1049        glb_error("Failure while copying input of glbChiNP");
1050        return -1;
1051      }
1052  }
1053 
1054  for(i=0;i<GLB_OSCP;i++) x[i]=glbGetOscParams(in,i);
1055  for(i=0;i<glb_num_of_exps;i++) x[i+GLB_OSCP]=glbGetDensityParams(in,i);
1056
1057  for(i=0;i<6+glb_num_of_exps;i++) fix_params[i]=x[i];
1058  for(i=0;i<n_free;i++) sp2[i+1]=x[index_tab[i]];
1059  if (setjmp(env)==1) 
1060    {
1061      okay_flag=1;
1062      return er1; 
1063    }
1064  glb_powell2(sp2,mat2,dim,TOLOSC,&it,&er1,MD_chi_NP);
1065  if(out!=NULL)
1066    {
1067      for(i=0;i<n_free;i++) 
1068        {
1069          if(index_tab[i]<GLB_OSCP)
1070            glbSetOscParams(out,sp2[i+1],index_tab[i]);
1071          else
1072            glbSetDensityParams(out,sp2[i+1],index_tab[i]-GLB_OSCP);
1073        }
1074      out=glbSetIteration(out,count);
1075    }
1076  glb_free_vec(sp2,1,dim);
1077  glb_free_mat(mat2,1,dim,1,dim); 
1078  return er1;
1079}
1080
1081double glbChiNP(const glb_params in, glb_params out, int exp)
1082{
1083  double res;
1084 
1085  if(in==NULL)
1086    {
1087      glb_error("Failure in glbChiNP: Input pointer must be non-NULL");
1088      return -1;
1089    }
1090 
1091  if(exp==GLB_ALL) res=internal_glbChiNP(in,out);
1092  else res=internal_glbSingleChiNP(in,out,exp);
1093  return res; 
1094}
1095
1096
1097/* Convenience wrappers for glbChiNP
1098 *
1099 * First the 1-d wrappers
1100 */
1101
1102static double glbChi1P(const glb_params in, 
1103                      glb_params out, int exp, int n)
1104{
1105  double res;
1106  int i,*b;
1107  int swit[37]; 
1108  int buff[37];
1109  b=CheckProjection();
1110  for(i=0;i<6+glb_num_of_exps;i++) buff[i]=b[i];
1111  for(i=0;i<6+glb_num_of_exps;i++) swit[i]=GLB_FREE;
1112  swit[n-1]=GLB_FIXED;
1113  SelectProjection(swit);
1114  if(in==NULL)
1115    {
1116      glb_error("Failure in glbChiNP: Input pointer must be non-NULL");
1117      return -1;
1118    }
1119 
1120  if(exp==GLB_ALL) res=internal_glbChiNP(in,out);
1121  else res=internal_glbSingleChiNP(in,out,exp);
1122  SelectProjection(buff);
1123  return res; 
1124}
1125
1126
1127double glbChiTheta(const glb_params in,glb_params out, int exp)
1128{
1129  double res;
1130  res=glbChi1P(in,out,exp,2);
1131  return res;
1132}
1133
1134double glbChiTheta23(const glb_params in,glb_params out, int exp)
1135{
1136  double res;
1137  res=glbChi1P(in,out,exp,3);
1138  return res;
1139}
1140
1141double glbChiDelta(const glb_params in,glb_params out, int exp)
1142{
1143  double res;
1144  res=glbChi1P(in,out,exp,4);
1145  return res;
1146}
1147
1148double glbChiDms(const glb_params in,glb_params out, int exp)
1149{
1150  double res;
1151  res=glbChi1P(in,out,exp,5);
1152  return res;
1153}
1154
1155double glbChiDm(const glb_params in,glb_params out, int exp)
1156{
1157  double res;
1158  res=glbChi1P(in,out,exp,6);
1159  return res;
1160}
1161
1162/* 2-d wrapper -- only ChiThetaDelta */
1163
1164double glbChiThetaDelta(const glb_params in,glb_params out, int exp)
1165{
1166  double res;
1167  int i,*b;
1168  int swit[37]; 
1169  int buff[37];
1170
1171  b=CheckProjection();
1172  for(i=0;i<6+glb_num_of_exps;i++) buff[i]=b[i];
1173  for(i=0;i<6+glb_num_of_exps;i++) swit[i]=GLB_FREE;
1174  swit[1]=GLB_FIXED;
1175  swit[3]=GLB_FIXED;
1176  SelectProjection(swit);
1177  if(in==NULL)
1178    {
1179      glb_error("Failure in glbChiNP: Input pointer must be non-NULL");
1180      return -1;
1181    }
1182 
1183  if(exp==GLB_ALL) res=internal_glbChiNP(in,out);
1184  else res=internal_glbSingleChiNP(in,out,exp);
1185  SelectProjection(b);
1186  return res;
1187}
1188
1189/* Wrapper for ChiAll */
1190
1191double glbChiAll(const glb_params in,glb_params out, int exp)
1192{
1193  double res;
1194  int i,*b;
1195  int swit[37]; 
1196  int buff[37];
1197  b=CheckProjection();
1198  for(i=0;i<6+glb_num_of_exps;i++) buff[i]=b[i];
1199  for(i=0;i<6+glb_num_of_exps;i++) swit[i]=GLB_FREE;
1200  SelectProjection(swit);
1201  if(in==NULL)
1202    {
1203      glb_error("Failure in glbChiNP: Input pointer must be non-NULL");
1204      return -1;
1205    }
1206 
1207  if(exp==GLB_ALL) res=internal_glbChiNP(in,out);
1208  else res=internal_glbSingleChiNP(in,out,exp);
1209  SelectProjection(b);
1210  return res;
1211}
1212
1213
1214/* New functions dealing with setting the projection flags */
1215
1216
1217int 
1218glbSetProjection(const glb_projection in)
1219{
1220  int i;
1221  int *buf;
1222  if(in==NULL) return -1;
1223  buf=glb_malloc(sizeof(int)*(GLB_OSCP+glb_num_of_exps));
1224
1225  for(i=0;i<GLB_OSCP;i++) buf[i]=glbGetProjectionFlag(in,i);
1226  for(i=0;i<glb_num_of_exps;i++) buf[i+GLB_OSCP]=
1227                                   glbGetDensityProjectionFlag(in,i);
1228  SelectProjection(buf);
1229  glb_free(buf);
1230  return 0;
1231
1232}
1233
1234int
1235glbGetProjection(glb_projection in)
1236{
1237  int i;
1238  if(in==NULL) return -1;
1239  for(i=0;i<GLB_OSCP;i++) in=glbSetProjectionFlag(in,para_tab[i],i);
1240  for(i=0;i<glb_num_of_exps;i++) 
1241    in=glbSetDensityProjectionFlag(in,para_tab[i+GLB_OSCP],i);
1242  return 0;
1243}
Note: See TracBrowser for help on using the repository browser.