GIRAFFE Pipeline Reference Manual

gimath_lm.c
1 /* $Id$
2  *
3  * This file is part of the GIRAFFE Pipeline
4  * Copyright (C) 2002-2006 European Southern Observatory
5  *
6  * This program is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * This program is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this program; if not, write to the Free Software
18  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
19  */
20 
21 /*
22  * $Author$
23  * $Date$
24  * $Revision$
25  * $Name$
26  */
27 
28 #ifdef HAVE_CONFIG_H
29 # include <config.h>
30 #endif
31 
32 #include <stdlib.h>
33 #include <stdio.h>
34 #include <math.h>
35 
36 #include <cxmemory.h>
37 #include <cxmessages.h>
38 
39 #include "gimacros.h"
40 #include "gimath.h"
41 #include "gimath_lm.h"
42 #include "gimessages.h"
43 
58 #define SWAP(a,b) {swap=(a);(a)=(b);(b)=swap;}
59 
60 
68 lmrq_model lmrq_models[] = {
69  {LMRQ_GAUSSUM, mrqgaussum, 4, 1, "gaussum", LINE_MODEL},
70  {LMRQ_XOPTMOD, mrqxoptmod, 7, 3, "xoptmod", XOPT_MODEL},
71  {LMRQ_XOPTMODGS, mrqxoptmodGS, 7, 3, "xoptmodGS", XOPT_MODEL},
72  {LMRQ_XOPTMOD2, mrqxoptmod2, 10, 3, "xoptmod2", XOPT_MODEL},
73  {LMRQ_PSFCOS, mrqpsfcos, 5, 1, "psfcos", LINE_MODEL},
74  {LMRQ_PSFEXP, mrqpsfexp, 5, 1, "psfexp", LINE_MODEL},
75  {LMRQ_YOPTMOD, mrqyoptmod, 7, 3, "yoptmod", YOPT_MODEL},
76  {LMRQ_YOPTMOD2, mrqyoptmod2, 10, 3, "yoptmod2", YOPT_MODEL},
77  {LMRQ_LOCYWARP, mrqlocywarp, 5, 4, "locywarp", LOCY_MODEL},
78  {LMRQ_PSFEXP2, mrqpsfexp2, 5, 1, "psfexp2", LINE_MODEL},
79  {LMRQ_TEST, mrqtest, 2, 1, "test", LINE_MODEL}
80 };
81 
82 cxint nr_lmrq_models = CX_N_ELEMENTS(lmrq_models);
83 
84 
101 static void
102 covariance_sort(cpl_matrix *covar, cxint ma, cxint ia[], cxint mfit)
103 {
104 
105  register cxint i, j, k;
106  register cxdouble swap;
107 
108  cxdouble *pd_covar = NULL;
109  cxint nr_covar;
110 
111  pd_covar = cpl_matrix_get_data(covar);
112  nr_covar = cpl_matrix_get_nrow(covar);
113 
114  for (i = mfit; i < ma; i++)
115  for (j = 0; j <= i; j++)
116  pd_covar[i * nr_covar + j] = pd_covar[j * nr_covar + i] = 0.0;
117 
118  k = mfit - 1;
119  for (j = (ma-1); j >= 0; j--) {
120  if (ia[j]) {
121  for (i = 0; i < ma; i++)
122  SWAP(pd_covar[i * nr_covar + k],pd_covar[i * nr_covar + j])
123 
124  for (i = 0;i < ma; i++)
125  SWAP(pd_covar[k * nr_covar + i],pd_covar[j * nr_covar + i])
126 
127  k--;
128  }
129  }
130 
131 } /* end covariance_sort() */
132 
133 #undef SWAP
134 
161 static double
162 mrqdydaweight(cxdouble x, cxdouble x0, cxdouble dx)
163 {
164  register cxdouble w;
165 
166  w = exp(-pow(fabs(x-x0),DW_DEGREE)/pow(dx,DW_DEGREE/DW_LOG001));
167 
168  if (isnan(w))
169  w = 1;
170 
171  return w;
172 }
173 
200 cxint
201 mrqnlfit(
202  cpl_matrix *x,
203  cpl_matrix *y,
204  cpl_matrix *sig,
205  cxint ndata,
206  cpl_matrix *a,
207  cxdouble r[],
208  cxint ia[],
209  cxint ma,
210  cpl_matrix *alpha,
211  cxdouble *chisq,
212  lmrq_params fit_params,
213  fitted_func funcs
214 ) {
215 
216  cxint itst,
217  n,
218  res;
219 
220  cxdouble alamda,
221  ochisq;
222 
223  cpl_matrix *beta = NULL;
224 
225  /*************************************************************************
226  PROCESSING
227  *************************************************************************/
228 
229  beta = cpl_matrix_new(ma,ma);
230 
231  alamda = -1.0;
232 
233  res = mymrqmin(x, y, sig, ndata, a, r, ia, ma, alpha, beta, chisq,
234  funcs, &alamda);
235 
236  if (res != 0) {
237  cpl_matrix_delete(beta); beta = NULL;
238  return res;
239  }
240 
241  itst=0;
242 
243  for (n = 1; n <= fit_params.imax; n++) {
244 
245  ochisq = *chisq;
246 
247  res = mymrqmin(x, y, sig, ndata, a, r, ia, ma, alpha, beta, chisq,
248  funcs, &alamda);
249 
250  if (res!=0) {
251  cpl_matrix_delete(beta); beta = NULL;
252  return res;
253  }
254 
255  if (*chisq > ochisq)
256  itst=0;
257  else if (fabs(ochisq-*chisq) < fit_params.dchsq)
258  itst++;
259 
260  if (itst > fit_params.tmax)
261  break;
262  }
263 
264  /* get covariance matrix */
265  alamda=0.0;
266 
267  res = mymrqmin(x, y, sig, ndata, a, r, ia, ma, alpha, beta, chisq,
268  funcs, &alamda);
269 
270  if (res != 0) {
271  cpl_matrix_delete(beta); beta = NULL;
272  return res;
273  }
274 
275  cpl_matrix_delete(beta); beta = NULL;
276 
277  return n;
278 
279 } /* end mrqnlfit() */
280 
337 cxint
338 mymrqmin(
339  cpl_matrix *x,
340  cpl_matrix *y,
341  cpl_matrix *sig,
342  cxint ndata,
343  cpl_matrix *a,
344  cxdouble r[],
345  cxint ia[],
346  cxint ma,
347  cpl_matrix *covar,
348  cpl_matrix *alpha,
349  cxdouble *chisq,
350  fitted_func funcs,
351  cxdouble *alamda
352 ) {
353 
354  register cxint gj, j, k, l, m;
355 
356  static cxdouble *pd_a, *pd_covar, *pd_alpha;
357  static cxint nr_covar, nr_alpha, nr_moneda, mfit;
358 
359  static cpl_matrix *matry, *mbeta, *mda, *moneda;
360  static cxdouble *atry, *beta, *da, *oneda, ochisq;
361 
362  /*************************************************************************
363  PROCESSING
364  *************************************************************************/
365 
366  pd_a = cpl_matrix_get_data(a);
367  pd_covar = cpl_matrix_get_data(covar);
368  pd_alpha = cpl_matrix_get_data(alpha);
369  nr_covar = cpl_matrix_get_nrow(covar);
370  nr_alpha = cpl_matrix_get_nrow(alpha);
371 
372  if (*alamda<0.0) {
373 
374  matry = cpl_matrix_new(ma,1);
375  atry = cpl_matrix_get_data(matry);
376 
377  mbeta = cpl_matrix_new(ma,1);
378  beta = cpl_matrix_get_data(mbeta);
379 
380  mda = cpl_matrix_new(ma,1);
381  da = cpl_matrix_get_data(mda);
382 
383  for (mfit = 0, j = 0; j < ma; j++)
384  if (ia[j])
385  mfit++;
386 
387  moneda = cpl_matrix_new(1,mfit);
388  oneda = cpl_matrix_get_data(moneda);
389 
390  *alamda = 0.001;
391 
392  gj = mymrqcof(x, y, sig, ndata, a, r, ia, ma, alpha, mbeta,
393  chisq, funcs);
394 
395  if (gj != 0) {
396  cpl_matrix_delete(moneda); moneda = NULL; oneda = NULL;
397  cpl_matrix_delete(mda); mda = NULL; da = NULL;
398  cpl_matrix_delete(mbeta); mbeta = NULL; beta = NULL;
399  cpl_matrix_delete(matry); matry = NULL; atry = NULL;
400  return gj;
401  }
402 
403  ochisq = (*chisq);
404 
405  for (j = 0; j < ma; j++)
406  atry[j] = pd_a[j];
407 
408  }
409 
410  nr_moneda = cpl_matrix_get_nrow(moneda);
411 
412  for (j = -1, l = 0; l < ma; l++) {
413  if (ia[l]) {
414  for (j++, k = -1, m = 0; m < ma; m++) {
415  if (ia[m]) {
416  k++;
417  pd_covar[j * nr_covar + k] = pd_alpha[j * nr_alpha + k];
418  }
419  }
420 
421  pd_covar[j * nr_covar + j] =
422  pd_alpha[j * nr_alpha + j] * (1.0 + (*alamda));
423 
424  oneda[j * nr_moneda + 0] = beta[j];
425  }
426  }
427 
428  gj = giraffe_gauss_jordan(covar, mfit, moneda, 1);
429 
430  if (gj != 0) {
431  cpl_matrix_delete(moneda); moneda = NULL; oneda = NULL;
432  cpl_matrix_delete(mda); mda = NULL; da = NULL;
433  cpl_matrix_delete(mbeta); mbeta = NULL; beta = NULL;
434  cpl_matrix_delete(matry); matry = NULL; atry = NULL;
435  return gj;
436  }
437 
438  for (j = 0; j < mfit; j++)
439  da[j] = oneda[j * nr_moneda + 0];
440 
441  if (*alamda == 0.0) {
442  covariance_sort(covar, ma, ia, mfit);
443  cpl_matrix_delete(moneda); moneda = NULL; oneda = NULL;
444  cpl_matrix_delete(mda); mda = NULL; da = NULL;
445  cpl_matrix_delete(mbeta); mbeta = NULL; beta = NULL;
446  cpl_matrix_delete(matry); matry = NULL; atry = NULL;
447  return 0;
448  }
449 
450  for (j = -1, l = 0; l < ma; l++)
451  if (ia[l])
452  atry[l] = pd_a[l] + da[++j];
453 
454  gj = mymrqcof(x, y, sig, ndata, matry, r, ia, ma, covar, mda,
455  chisq, funcs);
456 
457  if (gj != 0) {
458  cpl_matrix_delete(moneda); moneda = NULL; oneda = NULL;
459  cpl_matrix_delete(mda); mda = NULL; da = NULL;
460  cpl_matrix_delete(mbeta); mbeta = NULL; beta = NULL;
461  cpl_matrix_delete(matry); matry = NULL; atry = NULL;
462  return gj;
463  }
464 
465  if (*chisq < ochisq) {
466 
467  *alamda *= 0.1;
468  ochisq = *chisq;
469 
470  for (j = -1, l = 0; l < ma; l++) {
471  if (ia[l]) {
472  for (j++, k = -1, m = 0; m < ma; m++) {
473  if (ia[m]) {
474  k++;
475  pd_alpha[j * nr_alpha + k] =
476  pd_covar[j * nr_covar + k];
477  }
478  }
479 
480  beta[j] = da[j];
481  pd_a[l] = atry[l];
482  }
483  }
484 
485  } else {
486  *alamda *= 10.0;
487  *chisq = ochisq;
488  }
489 
490  return 0;
491 
492 } /* end mymrqmin() */
493 
521 cxint
522 mymrqcof(
523  cpl_matrix *x,
524  cpl_matrix *y,
525  cpl_matrix *sig,
526  cxint ndata,
527  cpl_matrix *a,
528  cxdouble r[],
529  cxint ia[],
530  cxint ma,
531  cpl_matrix *alpha,
532  cpl_matrix *beta,
533  cxdouble *chisq,
534  fitted_func funcs
535 ) {
536 
537  register cxint i, j, k, l, m, mfit = 0;
538 
539  cxdouble ymod, wt, sig2i, dy, *dyda;
540 
541  cxdouble *pd_x = NULL,
542  *pd_y = NULL,
543  *pd_sig = NULL,
544  *pd_a = NULL,
545  *pd_alpha = NULL,
546  *pd_beta = NULL;
547 
548  cxint nr_alpha, nc_x;
549 
550  /************************************************************************
551  PROCESSING
552  ************************************************************************/
553 
554  pd_x = cpl_matrix_get_data(x);
555  nc_x = cpl_matrix_get_ncol(x);
556  pd_y = cpl_matrix_get_data(y);
557  pd_sig = cpl_matrix_get_data(sig);
558  pd_a = cpl_matrix_get_data(a);
559  pd_alpha = cpl_matrix_get_data(alpha);
560  nr_alpha = cpl_matrix_get_nrow(alpha);
561  pd_beta = cpl_matrix_get_data(beta);
562 
563  for (j = 0; j < ma; j++) {
564  if (ia[j])
565  mfit++;
566  }
567 
568  for (j = 0; j < mfit; j++) {
569  for (k = 0; k <= j; k++)
570  pd_alpha[j * nr_alpha + k] = 0.0;
571 
572  pd_beta[j] = 0.0;
573  }
574 
575  *chisq = 0.0;
576 
577  dyda = (cxdouble *) cx_calloc(ma, sizeof(cxdouble));
578 
579  for (i = 0; i < ndata; i++) {
580 
581  (*funcs)(&(pd_x[i*nc_x]), pd_a, r, &ymod, dyda, ma);
582 
583  if (pd_sig[i]==0.0) {
584  continue;
585  }
586 
587  sig2i = 1.0 / (pd_sig[i] * pd_sig[i]);
588 
589  dy = pd_y[i] - ymod;
590 
591  for (j = -1, l = 0; l < ma; l++) {
592 
593  if (ia[l]) {
594  wt = dyda[l] * sig2i;
595  for (j++, k = -1, m = 0; m <= l; m++) {
596  if (ia[m]) {
597  ++k;
598  pd_alpha[j * nr_alpha + k] += (wt * dyda[m]);
599  }
600  }
601 
602  pd_beta[j] += (dy * wt);
603 
604  }
605  }
606 
607  *chisq += (dy * dy * sig2i);
608 
609  }
610 
611  for (j = 1; j < mfit; j++)
612  for (k = 0; k < j; k++)
613  pd_alpha[k * nr_alpha + j] = pd_alpha[j * nr_alpha + k];
614 
615 
616  cx_free(dyda);
617 
618  return 0;
619 
620 } /* end mymrqcof() */
621 
638 cxdouble
639 r_squared(cxdouble resSS, cpl_matrix *y, cxint n)
640 {
641  register cxint i;
642  register cxdouble Sy, Syy, SS;
643  cxdouble res, *pd_y = NULL;
644 
645  pd_y = cpl_matrix_get_data(y);
646 
647  if (n < 1)
648  return 0.0;
649 
650  for (i=0, Sy=0.0, Syy=0.0; i<n; i++) {
651  Sy += pd_y[i];
652  Syy += pd_y[i]*pd_y[i];
653  }
654 
655  SS = Syy - Sy*Sy/n;
656  res = resSS/SS;
657 
658  if (isnan(res))
659  return 0.0;
660 
661  if (res > 0.0)
662  res = sqrt(res);
663 
664  return res;
665 
666 } /* end r_squared() */
667 
699 void
700 mrqgaussum(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
701  cxdouble dyda[], cxint na)
702 {
703 
704  register cxint i, j;
705  register cxdouble fac,ex,amplitude,center,backg,width,xred;
706 
707  *y = 0.0;
708 
709  for (j = 0, i = 0; i < na; i += 4, j += 4) {
710  amplitude = a[i];
711  center = a[i + 1];
712  backg = a[i + 2];
713  width = a[i + 3];
714  xred = (x[0] - center) / width;
715  ex = exp(-xred * xred / 2.);
716  fac = amplitude * xred * ex;
717  *y += (amplitude * ex + backg);
718 
719  /* Check if derivatives expected */
720  if (dyda == NULL) continue;
721 
722  /* derivatives for each parameters */
723  dyda[j] = ex; /* d(y)/d(amplitude) */
724  dyda[j + 1] = fac / width; /* d(y)/d(center) */
725  dyda[j + 2] = 1.; /* d(y)/d(backg) */
726  dyda[j + 3] = (fac * xred) / width; /* d(y)/d(width) */
727  }
728 
729 } /* end mrqgaussum() */
730 
782 void
783 mrqxoptmod(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
784  cxdouble dyda[], cxint na)
785 {
786 
787  const cxchar *fctid = "mrqxoptmod";
788 
789  register cxdouble xccd, d, X;
790  register cxdouble lambda,xfibre,yfibre,pixsize,nx;
791  /* Optical model parameters */
792  register cxdouble fcoll,cfact;
793  /* Grating parameters */
794  register cxdouble gtheta,gorder,gspace;
795  register cxdouble yfibre2,tmp,tmp2,d2,X2,gspace2,sqtmp,costheta,sintheta;
796 
797  /* check for number of parameters */
798  if (na != 7) {
799  cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
800  return;
801  }
802 
803  *y = 0.0;
804  if (dyda != NULL) {
805  dyda[0] = dyda[1] = dyda[2] = dyda[3] =
806  dyda[4] = dyda[5] = dyda[6] = 0.0;
807  }
808 
809  lambda = x[LMI_WLEN]; /* wavelength [mm] */
810  xfibre = x[LMI_XFIB]; /* Y fibre [mm] */
811  yfibre = x[LMI_YFIB]; /* Y fibre [mm] */
812 
813  nx = a[LMP_NX]; /* CCD size in X [pixels] */
814  pixsize = a[LMP_PXSIZ]; /* CCD pixel size [mm] */
815  fcoll = a[LMP_FCOLL]; /* collimator focal length [mm] */
816  cfact = a[LMP_CFACT]; /* camera magnification factor */
817  gtheta = a[LMP_THETA]; /* grating angle [radian] */
818  gorder = a[LMP_ORDER]; /* grating diffraction order */
819  gspace = a[LMP_SPACE]; /* grating groove spacing [mm] */
820 
821  yfibre2 = yfibre * yfibre;
822  gspace2 = gspace * gspace;
823  costheta = cos(gtheta);
824  sintheta = sin(gtheta);
825  d2 = xfibre * xfibre + yfibre2 + (fcoll * fcoll);
826  d = sqrt(d2);
827  X = (-lambda*gorder/gspace) + (xfibre*costheta/d) + (fcoll*sintheta/d);
828  X2 = X * X;
829  sqtmp = sqrt(1.0 - yfibre2/d2 - X2);
830  tmp = -sintheta*X + costheta*sqtmp;
831  tmp2 = tmp * tmp;
832  xccd = (cfact * fcoll * (X*costheta + sintheta*sqtmp))/tmp;
833 
834  /* takes care of model direction */
835  if (nx < 0.0)
836  *y = (xccd / pixsize - 0.5*nx);
837  else
838  *y = (-xccd / pixsize + 0.5*nx);
839 
840  /* Check if derivatives expected */
841  if (dyda == NULL)
842  return;
843 
844  /* derivatives for each parameters */
845  dyda[LMP_NX] = 0.5; /* d(y)/d(nx) */
846  dyda[LMP_PXSIZ] = 0.0; /* d(y)/d(pixsize) */
847 
848  dyda[LMP_FCOLL] = cfact*(costheta*X+sintheta*sqtmp)/tmp +
849  cfact*fcoll*(costheta*(-X*fcoll/d2+sintheta/d -
850  gorder*lambda*fcoll/(d2*gspace)) +
851  0.5*sintheta*(-2.0*X*(-X*fcoll/d2+sintheta/d -
852  gorder*lambda*fcoll/(d2*gspace))+2.0*yfibre2*fcoll/(d2*d2))/sqtmp)/tmp -
853  cfact*fcoll*(costheta*X+sintheta*sqtmp)*(-sintheta*(-X*fcoll/d2 +
854  sintheta/d-gorder*lambda*fcoll/(d2*gspace)) +
855  0.5*costheta*(-2.0*X*(-X*fcoll/d2+sintheta/d -
856  gorder*lambda*fcoll/(d2*gspace))+2.0*yfibre2*fcoll/(d2*d2))/sqtmp)/tmp2;
857  dyda[LMP_FCOLL] /= pixsize; /* d(y)/d(fcoll) */
858 
859  dyda[LMP_CFACT] = (xccd/cfact)/pixsize; /* d(y)/d(cfact) */
860 
861  dyda[LMP_THETA] = cfact*fcoll*((-xfibre*sintheta/d+fcoll*costheta/d)*costheta -
862  sintheta*X-sintheta*X*(-xfibre*sintheta/d+fcoll*costheta/d)/sqtmp +
863  costheta*sqtmp)/tmp -
864  cfact*fcoll*(costheta*X+sintheta*sqtmp)*(-(-xfibre*sintheta/d +
865  fcoll*costheta/d)*sintheta-costheta*X -
866  costheta*X*(-xfibre*sintheta/d+fcoll*costheta/d)/sqtmp -
867  sintheta*sqtmp)/tmp2;
868  dyda[LMP_THETA] /= pixsize; /* d(y)/d(gtheta) */
869 
870  dyda[LMP_ORDER] = 0.0; /* d(y)/d(gorder) */
871  dyda[LMP_SPACE] = cfact*fcoll*(lambda*gorder*costheta/gspace2-sintheta*X*lambda*gorder/(sqtmp*gspace2))/tmp -
872  cfact*fcoll*(X*costheta+sintheta*sqtmp) *
873  (-lambda*gorder*sintheta/gspace2-costheta*X*lambda*gorder/(sqtmp*gspace2))/tmp2;
874  dyda[LMP_SPACE] /= pixsize; /* d(y)/d(gspace) */
875 
876  if (nx > 0.0) {
877  dyda[LMP_NX] = -dyda[LMP_NX];
878  dyda[LMP_PXSIZ] = -dyda[LMP_PXSIZ];
879  dyda[LMP_FCOLL] = -dyda[LMP_FCOLL];
880  dyda[LMP_CFACT] = -dyda[LMP_CFACT];
881  dyda[LMP_THETA] = -dyda[LMP_THETA];
882  dyda[LMP_ORDER] = -dyda[LMP_ORDER];
883  dyda[LMP_SPACE] = -dyda[LMP_SPACE];
884  }
885 
886  if (r != NULL) {
887  register cxint k;
888 
889  k = LMP_FCOLL << 1;
890  if (r[k+1] > 0) {
891  dyda[LMP_FCOLL] *= mrqdydaweight(a[LMP_FCOLL],r[k],r[k+1]);
892  }
893  k = LMP_CFACT << 1;
894  if (r[k+1] > 0) {
895  dyda[LMP_CFACT] *= mrqdydaweight(a[LMP_CFACT],r[k],r[k+1]);
896  }
897  k = LMP_THETA << 1;
898  if (r[k+1] > 0) {
899  dyda[LMP_THETA] *= mrqdydaweight(a[LMP_THETA],r[k],r[k+1]);
900  }
901  k = LMP_SPACE << 1;
902  if (r[k+1] > 0) {
903  dyda[LMP_SPACE] *= mrqdydaweight(a[LMP_SPACE],r[k],r[k+1]);
904  }
905  }
906 
907 } /* end mrqxoptmod() */
908 
970 void
971 mrqxoptmod2(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
972  cxdouble dyda[], cxint na)
973 {
974 
975  const cxchar *fctid = "mrqxoptmod2";
976 
977  register cxdouble lambda,xfibre,yfibre,pixsize,nx;
978  /* Optical model parameters */
979  register cxdouble fcoll,cfact;
980  /* Grating parameters */
981  register cxdouble gtheta,gorder,gspace;
982  /* Slit position parameters */
983  cxdouble slitdx,slitdy,slitphi;
984 
985  register cxdouble t1,t10,t104,t107,t11,t113,t119,t12,t120,t121,t124,t136,
986  t137,t138,t14,t143,t148,t16,t161,t162,t166,t168,t17,t173,
987  t18,t19,t191,t195,t196,t2,t20,t201,t21,t210,t23,t24,t26,
988  t27,t28,t3,t30,t32,t33,t34,t35,t36,t37,t38,t39,t4,t40,t44,
989  t49,t52,t58,t60,t61,t62,t64,t68,t75,t76,t78,t80,t9,t91,t93;
990 
991  /* check for number of parameters */
992  if (na != 10) {
993  cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
994  return;
995  }
996 
997  *y = 0.0;
998  if (dyda != NULL) {
999  dyda[0] = dyda[1] = dyda[2] = dyda[3] =
1000  dyda[4] = dyda[5] = dyda[6] =
1001  dyda[7] = dyda[8] = dyda[9] = 0.0;
1002  }
1003 
1004  lambda = x[LMI_WLEN]; /* wavelength [mm] */
1005  xfibre = x[LMI_XFIB]; /* Y fibre [mm] */
1006  yfibre = x[LMI_YFIB]; /* Y fibre [mm] */
1007 
1008  nx = a[LMP_NX]; /* CCD size in X [pixels] */
1009  pixsize = a[LMP_PXSIZ]; /* CCD pixel size [mm] */
1010  fcoll = a[LMP_FCOLL]; /* collimator focal length [mm] */
1011  cfact = a[LMP_CFACT]; /* camera magnification factor */
1012  gtheta = a[LMP_THETA]; /* grating angle [radian] */
1013  gorder = a[LMP_ORDER]; /* grating diffraction order */
1014  gspace = a[LMP_SPACE]; /* grating groove spacing [mm] */
1015  slitdx = a[LMP_SOFFX]; /* slit position x offset [mm] */
1016  slitdy = a[LMP_SOFFY]; /* slit position y offset [mm] */
1017  slitphi = a[LMP_SPHI]; /* slit position angle [radian] */
1018 
1019  t1 = cfact*fcoll;
1020  t2 = cos(gtheta);
1021  t3 = lambda*gorder;
1022  t4 = 1.0/gspace;
1023  t9 = xfibre*(1.0+slitphi*yfibre)+slitdx;
1024  t10 = t9*t2;
1025  t11 = t9*t9;
1026  t12 = slitphi*slitphi;
1027  t14 = sqrt(1.0-t12);
1028  t16 = yfibre*t14+slitdy;
1029  t17 = t16*t16;
1030  t18 = fcoll*fcoll;
1031  t19 = t11+t17+t18;
1032  t20 = sqrt(t19);
1033  t21 = 1.0/t20;
1034  t23 = sin(gtheta);
1035  t24 = fcoll*t23;
1036  t26 = -t3*t4+t10*t21+t24*t21;
1037  t27 = t2*t26;
1038  t28 = 1.0/t19;
1039  t30 = t26*t26;
1040  t32 = sqrt(1.0-t17*t28-t30);
1041  t33 = t23*t32;
1042  t34 = t27+t33;
1043  t35 = t23*t26;
1044  t36 = t2*t32;
1045  t37 = -t35+t36;
1046  t38 = 1.0/t37;
1047  t39 = t34*t38;
1048  t40 = 1.0/pixsize;
1049  t44 = pixsize*pixsize;
1050  t49 = t38*t40;
1051  t52 = 1.0/t20/t19;
1052  t58 = -t10*t52*fcoll+t23*t21-t18*t23*t52;
1053  t60 = 1.0/t32;
1054  t61 = t23*t60;
1055  t62 = t19*t19;
1056  t64 = t17/t62;
1057  t68 = 2.0*t64*fcoll-2.0*t26*t58;
1058  t75 = t1*t34;
1059  t76 = t37*t37;
1060  t78 = 1.0/t76*t40;
1061  t80 = t2*t60;
1062  t91 = -t9*t23*t21+fcoll*t2*t21;
1063  t93 = t26*t91;
1064  t104 = t2*lambda;
1065  t107 = t26*lambda*t4;
1066  t113 = t23*lambda;
1067  t119 = gspace*gspace;
1068  t120 = 1.0/t119;
1069  t121 = gorder*t120;
1070  t124 = t3*t120;
1071  t136 = t2*t21;
1072  t137 = 2.0*t9;
1073  t138 = t52*t137;
1074  t143 = t136-t10*t138/2.0-t24*t138/2.0;
1075  t148 = t64*t137-2.0*t26*t143;
1076  t161 = 2.0*t16;
1077  t162 = t52*t161;
1078  t166 = -t10*t162/2.0-t24*t162/2.0;
1079  t168 = t16*t28;
1080  t173 = -2.0*t168+t64*t161-2.0*t26*t166;
1081  t191 = 1.0/t14;
1082  t195 = 2.0*t9*xfibre*yfibre-2.0*t16*yfibre*t191*slitphi;
1083  t196 = t52*t195;
1084  t201 = xfibre*yfibre*t136-t10*t196/2.0-t24*t196/2.0;
1085  t210 = 2.0*t168*yfibre*t191*slitphi+t64*t195-2.0*t26*t201;
1086 
1087  /* takes care of model direction */
1088  if (nx < 0.0)
1089  *y = t1*t39*t40-0.5*nx;
1090  else
1091  *y = -t1*t39*t40+0.5*nx;
1092 
1093  /* Check if derivatives expected */
1094  if (dyda == NULL)
1095  return;
1096 
1097  /* derivatives for each parameters */
1098  dyda[LMP_NX] = 0.5; /* d(y)/d(nx) */
1099  dyda[LMP_PXSIZ] = -t1*t39/t44; /* d(y)/d(pixsize) */
1100  dyda[LMP_FCOLL] = /* d(y)/d(fcoll) */
1101  cfact*t34*t49+t1*(t2*t58+t61*t68/2.0)*t38*t40 -
1102  t75*t78*(-t23*t58+t80*t68/2.0);
1103  dyda[LMP_CFACT] = /* d(y)/d(cfact) */
1104  fcoll*t34*t49;
1105  dyda[LMP_THETA] = /* d(y)/d(gtheta) */
1106  t1*(-t35+t2*t91+t36-t61*t93)*t38*t40 -
1107  t75*t78*(-t27-t23*t91-t33-t80*t93);
1108  dyda[LMP_ORDER] = /* d(y)/d(gorder) */
1109  t1*(-t104*t4+t61*t107)*t38*t40-t75*t78*(t113*t4+t80*t107);
1110  dyda[LMP_SPACE] = /* d(y)/d(gspace) */
1111  t1*(t104*t121-t61*t26*t124)*t38*t40 -
1112  t75*t78*(-t113*t121-t80*t26*t124);
1113  dyda[LMP_SOFFX] = /* d(y)/d(slitdx) */
1114  t1*(t2*t143+t61*t148/2.0)*t38*t40 -
1115  t75*t78*(-t23*t143+t80*t148/2.0);
1116  dyda[LMP_SOFFY] = /* d(y)/d(slitdy) */
1117  t1*(t2*t166+t61*t173/2.0)*t38*t40 -
1118  t75*t78*(-t23*t166+t80*t173/2.0);
1119  dyda[LMP_SPHI] = /* d(y)/d(slitphi) */
1120  t1*(t2*t201+t61*t210/2.0)*t38*t40 -
1121  t75*t78*(-t23*t201+t80*t210/2.0);
1122 
1123  if (nx > 0.0) {
1124  dyda[LMP_NX] = -dyda[LMP_NX];
1125  dyda[LMP_PXSIZ] = -dyda[LMP_PXSIZ];
1126  dyda[LMP_FCOLL] = -dyda[LMP_FCOLL];
1127  dyda[LMP_CFACT] = -dyda[LMP_CFACT];
1128  dyda[LMP_THETA] = -dyda[LMP_THETA];
1129  dyda[LMP_ORDER] = -dyda[LMP_ORDER];
1130  dyda[LMP_SPACE] = -dyda[LMP_SPACE];
1131  dyda[LMP_SOFFX] = -dyda[LMP_SOFFX];
1132  dyda[LMP_SOFFY] = -dyda[LMP_SOFFY];
1133  dyda[LMP_SPHI] = -dyda[LMP_SPHI];
1134  }
1135 
1136  if (r != NULL) {
1137  register cxint k;
1138 
1139  k = LMP_PXSIZ << 1;
1140  if (r[k+1] > 0) {
1141  dyda[LMP_PXSIZ] *= mrqdydaweight(a[LMP_PXSIZ],r[k],r[k+1]);
1142  }
1143  k = LMP_FCOLL << 1;
1144  if (r[k+1] > 0) {
1145  dyda[LMP_FCOLL] *= mrqdydaweight(a[LMP_FCOLL],r[k],r[k+1]);
1146  }
1147  k = LMP_CFACT << 1;
1148  if (r[k+1] > 0) {
1149  dyda[LMP_CFACT] *= mrqdydaweight(a[LMP_CFACT],r[k],r[k+1]);
1150  }
1151  k = LMP_THETA << 1;
1152  if (r[k+1] > 0) {
1153  dyda[LMP_THETA] *= mrqdydaweight(a[LMP_THETA],r[k],r[k+1]);
1154  }
1155  k = LMP_ORDER << 1;
1156  if (r[k+1] > 0) {
1157  dyda[LMP_ORDER] *= mrqdydaweight(a[LMP_ORDER],r[k],r[k+1]);
1158  }
1159  k = LMP_SPACE << 1;
1160  if (r[k+1] > 0) {
1161  dyda[LMP_SPACE] *= mrqdydaweight(a[LMP_SPACE],r[k],r[k+1]);
1162  }
1163  k = LMP_SOFFX << 1;
1164  if (r[k+1] > 0) {
1165  dyda[LMP_SOFFX] *= mrqdydaweight(a[LMP_SOFFX],r[k],r[k+1]);
1166  }
1167  k = LMP_SOFFY << 1;
1168  if (r[k+1] > 0) {
1169  dyda[LMP_SOFFY] *= mrqdydaweight(a[LMP_SOFFY],r[k],r[k+1]);
1170  }
1171  k = LMP_SPHI << 1;
1172  if (r[k+1] > 0) {
1173  dyda[LMP_SPHI] *= mrqdydaweight(a[LMP_SPHI],r[k],r[k+1]);
1174  }
1175  }
1176 
1177 } /* end mrqxoptmod2() */
1178 
1243 void
1244 mrqyoptmod(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
1245  cxdouble dyda[], cxint na)
1246 {
1247 
1248  const cxchar *fctid = "mrqyoptmod";
1249 
1250  register cxdouble lambda,xfibre,yfibre,pixsize,ny;
1251  /* Optical model parameters */
1252  register cxdouble fcoll,cfact;
1253  /* Grating parameters */
1254  register cxdouble gtheta,gorder,gspace;
1255 
1256  cxdouble t10,t12,t13,t15,t18,t2,t22,t24,t26,t27,t28,t29,t3,t30,t33,
1257  t4,t41,t45,t47,t5,t53,t56,t57,t6,t7,t76,t8,t9,t93,t94;
1258 
1259  /* check for number of parameters */
1260  if (na != 7) {
1261  cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
1262  return;
1263  }
1264 
1265  *y = 0.0;
1266  if (dyda != NULL) {
1267  dyda[LMP_NX] = dyda[LMP_PXSIZ] = dyda[LMP_FCOLL] = dyda[LMP_CFACT] =
1268  dyda[LMP_THETA] = dyda[LMP_ORDER] = dyda[LMP_SPACE] = 0.0;
1269  }
1270 
1271  lambda = x[LMI_WLEN];
1272  xfibre = x[LMI_XFIB];
1273  yfibre = x[LMI_YFIB];
1274 
1275  ny = a[LMP_NY];
1276  pixsize = a[LMP_PXSIZ];
1277  fcoll = a[LMP_FCOLL];
1278  cfact = a[LMP_CFACT];
1279  gtheta = a[LMP_THETA];
1280  gorder = a[LMP_ORDER];
1281  gspace = a[LMP_SPACE];
1282 
1283  t2 = cfact*fcoll*yfibre;
1284  t3 = xfibre*xfibre;
1285  t4 = yfibre*yfibre;
1286  t5 = fcoll*fcoll;
1287  t6 = t3+t4+t5;
1288  t7 = sqrt(t6);
1289  t8 = 1.0/t7;
1290  t9 = lambda*gorder;
1291  t10 = 1.0/gspace;
1292  t12 = cos(gtheta);
1293  t13 = xfibre*t12;
1294  t15 = sin(gtheta);
1295  t18 = -t9*t10+t13*t8+fcoll*t15*t8;
1296  t22 = t18*t18;
1297  t24 = sqrt(1.0-t4/t6-t22);
1298  t26 = -t18*t15+t12*t24;
1299  t27 = 1.0/t26;
1300  t28 = t8*t27;
1301  t29 = 1.0/pixsize;
1302  t30 = t28*t29;
1303  t33 = pixsize*pixsize;
1304  t41 = 1.0/t7/t6;
1305  t45 = t26*t26;
1306  t47 = t8/t45;
1307  t53 = -t13*t41*fcoll+t15*t8-t5*t15*t41;
1308  t56 = t12/t24;
1309  t57 = t6*t6;
1310  t76 = -xfibre*t15*t8+fcoll*t12*t8;
1311  t93 = gspace*gspace;
1312  t94 = 1.0/t93;
1313 
1314  *y = -t2*t30+0.5*ny;
1315 
1316  /* Check if derivatives expected */
1317  if (dyda == NULL) return;
1318 
1319  /* derivatives for each parameters */
1320  dyda[LMP_NY] = 0.5; /* d(y)/d(ny) */
1321 
1322  dyda[LMP_PXSIZ] = t2*t28/t33;
1323  dyda[LMP_FCOLL] = -cfact*yfibre*t30+cfact*t5*yfibre*t41*t27*t29+t2*t47*t29*
1324  (-t53*t15+t56*(2.0*t4/t57*fcoll-2.0*t18*t53)/2.0);
1325  dyda[LMP_CFACT] = -fcoll*yfibre*t30;
1326  dyda[LMP_THETA] = t2*t47*t29*(-t76*t15-t18*t12-t15*t24-t56*t18*t76);
1327  dyda[LMP_ORDER] = t2*t47*t29*(lambda*t10*t15+t56*t18*lambda*t10);
1328  dyda[LMP_SPACE] = t2*t47*t29*(-t9*t94*t15-t56*t18*t9*t94);
1329 
1330 } /* end mrqyoptmod() */
1331 
1404 void
1405 mrqyoptmod2(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
1406  cxdouble dyda[], cxint na)
1407 {
1408 
1409  const cxchar *fctid = "mrqyoptmod2";
1410 
1411  register cxdouble lambda,xfibre,yfibre,pixsize,ny;
1412  /* Optical model parameters */
1413  register cxdouble fcoll,cfact;
1414  /* Grating parameters */
1415  register cxdouble gtheta,gorder,gspace;
1416  /* Slit position parameters */
1417  cxdouble slitdx,slitdy,slitphi;
1418 
1419  double t1,t102,t103,t11,t112,t117,t118,t12,t123,t13,t136,t14,t141,t145,
1420  t147,t15,t159,t16,t160,t17,t172,t179,t18,t184,t19,t2,t21,t22,t24,
1421  t25,t27,t29,t31,t33,t35,t36,t37,t38,t39,t4,t42,t50,t51,t54,t56,t6,
1422  t62,t65,t66,t68,t7,t85;
1423 
1424  /* check for number of parameters */
1425  if (na != 10) {
1426  cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
1427  return;
1428  }
1429 
1430  *y = 0.0;
1431  if (dyda != NULL) {
1432  dyda[LMP_NY] = dyda[LMP_PXSIZ] = dyda[LMP_FCOLL] = dyda[LMP_CFACT] =
1433  dyda[LMP_THETA] = dyda[LMP_ORDER] = dyda[LMP_SPACE] =
1434  dyda[LMP_SOFFX] = dyda[LMP_SOFFY] = dyda[LMP_SPHI] = 0.0;
1435  }
1436 
1437  lambda = x[LMI_WLEN];
1438  xfibre = x[LMI_XFIB];
1439  yfibre = x[LMI_YFIB];
1440 
1441  ny = a[LMP_NY];
1442  pixsize = a[LMP_PXSIZ];
1443  fcoll = a[LMP_FCOLL];
1444  cfact = a[LMP_CFACT];
1445  gtheta = a[LMP_THETA];
1446  gorder = a[LMP_ORDER];
1447  gspace = a[LMP_SPACE];
1448  slitdx = a[LMP_SOFFX];
1449  slitdy = a[LMP_SOFFY];
1450  slitphi = a[LMP_SPHI];
1451 
1452  t1 = cfact*fcoll;
1453  t2 = slitphi*slitphi;
1454  t4 = sqrt(1.0-t2);
1455  t6 = yfibre*t4+slitdy;
1456  t7 = t1*t6;
1457  t11 = xfibre*(1.0+slitphi*yfibre)+slitdx;
1458  t12 = t11*t11;
1459  t13 = t6*t6;
1460  t14 = fcoll*fcoll;
1461  t15 = t12+t13+t14;
1462  t16 = sqrt(t15);
1463  t17 = 1/t16;
1464  t18 = lambda*gorder;
1465  t19 = 1/gspace;
1466  t21 = cos(gtheta);
1467  t22 = t11*t21;
1468  t24 = sin(gtheta);
1469  t25 = fcoll*t24;
1470  t27 = -t18*t19+t22*t17+t25*t17;
1471  t29 = 1/t15;
1472  t31 = t27*t27;
1473  t33 = sqrt(1.0-t13*t29-t31);
1474  t35 = -t27*t24+t21*t33;
1475  t36 = 1/t35;
1476  t37 = t17*t36;
1477  t38 = 1/pixsize;
1478  t39 = t37*t38;
1479  t42 = pixsize*pixsize;
1480  t50 = 1/t16/t15;
1481  t51 = t50*t36;
1482  t54 = t35*t35;
1483  t56 = t17/t54;
1484  t62 = -t22*t50*fcoll+t24*t17-t14*t24*t50;
1485  t65 = t21/t33;
1486  t66 = t15*t15;
1487  t68 = t13/t66;
1488  t85 = -t11*t24*t17+fcoll*t21*t17;
1489  t102 = gspace*gspace;
1490  t103 = 1/t102;
1491  t112 = 2.0*t11;
1492  t117 = t21*t17;
1493  t118 = t50*t112;
1494  t123 = t117-t22*t118/2.0-t25*t118/2.0;
1495  t136 = 2.0*t6;
1496  t141 = t50*t136;
1497  t145 = -t22*t141/2.0-t25*t141/2.0;
1498  t147 = t6*t29;
1499  t159 = 1/t4;
1500  t160 = yfibre*t159;
1501  t172 = 2.0*t11*xfibre*yfibre-2.0*t6*yfibre*t159*slitphi;
1502  t179 = t50*t172;
1503  t184 = xfibre*yfibre*t117-t22*t179/2.0-t25*t179/2.0;
1504 
1505  *y = -t7*t39+0.5*ny;
1506 
1507  /* Check if derivatives expected */
1508  if (dyda == NULL) return;
1509 
1510  /* derivatives for each parameters */
1511  dyda[LMP_NY] = 0.5; /* d(y)/d(ny) */
1512  dyda[LMP_PXSIZ] = t7*t37/t42; /* d(y)/d(pixsize) */
1513  dyda[LMP_FCOLL] = /* d(y)/d(fcoll) */
1514  -cfact*t6*t39+cfact*t14*t6*t51*t38+
1515  t7*t56*t38*(-t62*t24+t65*(2.0*t68*fcoll-2.0*t27*t62)/2.0);
1516  dyda[LMP_CFACT] = /* d(y)/d(cfact) */
1517  -fcoll*t6*t39;
1518  dyda[LMP_THETA] = /* d(y)/d(gtheta) */
1519  t7*t56*t38*(-t85*t24-t27*t21-t24*t33-t65*t27*t85);
1520  dyda[LMP_ORDER] = /* d(y)/d(gorder) */
1521  t7*t56*t38*(lambda*t19*t24+t65*t27*lambda*t19);
1522  dyda[LMP_SPACE] = /* d(y)/d(gspace) */
1523  t7*t56*t38*(-t18*t103*t24-t65*t27*t18*t103);
1524  dyda[LMP_SOFFX] = /* d(y)/d(slitdx) */
1525  t7*t51*t38*t112/2.0+t7*t56*t38*(-t123*t24+t65*
1526  (t68*t112-2.0*t27*t123)/2.0);
1527  dyda[LMP_SOFFY] = /* d(y)/d(slitdy) */
1528  -t1*t39+t7*t51*t38*t136/2.0+t7*t56*t38*(-t145*t24+t65*
1529  (-2.0*t147+t68*t136-2.0*t27*t145)/2.0);
1530  dyda[LMP_SPHI] = /* d(y)/d(slitphi) */
1531  t1*t160*slitphi*t17*t36*t38+t7*t51*t38*t172/2.0+
1532  t7*t56*t38*(-t184*t24+t65*(2.0*t147*t160*slitphi+t68*t172-
1533  2.0*t27*t184)/2.0);
1534 
1535 } /* end mrqyoptmod2() */
1536 
1562 void
1563 mrqpsfcos(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
1564  cxdouble dyda[], cxint na)
1565 {
1566 
1567  const cxchar *fctid = "mrqpsfcos";
1568 
1569  cxdouble t1,t10,t13,t14,t15,t16,t2,t26,t3,t4,t5,t6,t7,t8,t9;
1570 
1571  cxdouble amplitude = a[LMP_AMPL];
1572  cxdouble center = a[LMP_CENT];
1573  cxdouble background = a[LMP_BACK];
1574  cxdouble width1 = a[LMP_WID1];
1575  cxdouble width2 = a[LMP_WID2];
1576 
1577  /* check for number of parameters */
1578  if (na != 5) {
1579  cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
1580  return;
1581  }
1582 
1583  *y = 0.0;
1584  if (dyda != NULL) {
1585  dyda[LMP_AMPL] = dyda[LMP_CENT] = dyda[LMP_BACK] = dyda[LMP_WID1] =
1586  dyda[LMP_WID2] = 0.0;
1587  }
1588 
1589  t1 = x[0]-center;
1590  t2 = fabs(t1);
1591  t3 = 1.0/width2;
1592  t4 = t2*t3;
1593  t5 = pow(t4,width1);
1594  t6 = CX_PI*t5;
1595  t7 = cos(t6);
1596  t8 = 1.0+t7;
1597  t9 = t8*t8;
1598  t10 = t9*t8;
1599  t13 = amplitude*t9;
1600  t14 = sin(t6);
1601  t15 = t13*t14;
1602  t16 = log(t4);
1603  t26 = (t1>0.0)? 1.0:-1.0;
1604 
1605  if (t2 > width2) {
1606  *y = background;
1607 
1608  /* Check if derivatives expected */
1609  if (dyda == NULL) return;
1610 
1611  dyda[LMP_AMPL] = dyda[LMP_CENT] = dyda[LMP_BACK] = dyda[LMP_WID1] = 0.0;
1612  dyda[LMP_WID2] = 1.0;
1613  } else {
1614  *y = amplitude*t10/8.0+background; /* Function value */
1615 
1616  /* Check if derivatives expected */
1617  if (dyda == NULL)
1618  return;
1619 
1620  dyda[LMP_AMPL] = t10/8.0; /* d(y)/d(amplitude) */
1621  /* d(y)/d(center) */
1622  dyda[LMP_CENT] = 3.0/8.0*t13*t14*CX_PI*t5*width1*t26/t2;
1623  dyda[LMP_BACK] = 1.0; /* d(y)/d(background) */
1624  dyda[LMP_WID1] = -3.0/8.0*t15*t6*t16; /* d(y)/d(width1) */
1625  dyda[LMP_WID2] = 3.0/8.0*t15*t6*width1*t3; /* d(y)/d(width2) */
1626  }
1627 
1628  return;
1629 
1630 } /* end mrqpsfcos() */
1631 
1657 void
1658 mrqpsfexp(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
1659  cxdouble dyda[], cxint na)
1660 {
1661 
1662  const cxchar *fctid = "mrqpsfexp";
1663 
1664  cxdouble t1,t2,t3,t4,t6,t8,t10,t15,t18;
1665 
1666  cxdouble amplitude = a[LMP_AMPL];
1667  cxdouble center = a[LMP_CENT];
1668  cxdouble background = a[LMP_BACK];
1669  cxdouble width1 = a[LMP_WID1];
1670  cxdouble width2 = a[LMP_WID2];
1671 
1672  /* check for number of parameters */
1673  if (na != 5) {
1674  cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
1675  return;
1676  }
1677 
1678  *y = 0.0;
1679  if (dyda != NULL) {
1680  dyda[LMP_AMPL] = dyda[LMP_CENT] = dyda[LMP_BACK] = dyda[LMP_WID1] =
1681  dyda[LMP_WID2] = 0.0;
1682  }
1683 
1684  t1 = x[0]-center;
1685 
1686  if (t1 > 0.0) {
1687  t2 = t1;
1688  t10 = 1.0;
1689  } else {
1690  t2 = -t1;
1691  t10 = -1.0;
1692  }
1693 
1694  t3 = pow(t2,width2);
1695  t4 = 1.0/width1;
1696  t6 = exp(-t3*t4);
1697  t8 = amplitude*t3;
1698  t15 = width1*width1;
1699  t18 = log(t2);
1700 
1701  *y = amplitude*t6+background;
1702 
1703  /* Check if derivatives expected */
1704  if (dyda == NULL)
1705  return;
1706 
1707  dyda[LMP_AMPL] = t6; /* d(y)/d(amplitude) */
1708  dyda[LMP_CENT] = t8*width2*t10/t2*t4*t6; /* d(y)/d(center) */
1709 
1710  if (isnan(dyda[LMP_CENT]))
1711  dyda[LMP_CENT] = 0.;
1712 
1713  dyda[LMP_BACK] = 1.0; /* d(y)/d(background) */
1714  dyda[LMP_WID1] = t8/t15*t6; /* d(y)/d(width1) */
1715  dyda[LMP_WID2] = -t8*t18*t4*t6; /* d(y)/d(width2) */
1716 
1717  if (isnan(dyda[LMP_WID2]))
1718  dyda[LMP_WID2] = 0.;
1719 
1720  if (r != NULL) {
1721  register cxint k;
1722 
1723  k = LMP_AMPL << 1;
1724  if (r[k+1] > 0) {
1725  dyda[LMP_AMPL] *= mrqdydaweight(a[LMP_AMPL],r[k],r[k+1]);
1726  }
1727  k = LMP_CENT << 1;
1728  if (r[k+1] > 0) {
1729  dyda[LMP_CENT] *= mrqdydaweight(a[LMP_CENT],r[k],r[k+1]);
1730  }
1731  k = LMP_WID1 << 1;
1732  if (r[k+1] > 0) {
1733  dyda[LMP_WID1] *= mrqdydaweight(a[LMP_WID1],r[k],r[k+1]);
1734  }
1735  k = LMP_WID2 << 1;
1736  if (r[k+1] > 0) {
1737  dyda[LMP_WID2] *= mrqdydaweight(a[LMP_WID2],r[k],r[k+1]);
1738  }
1739  }
1740 
1741  return;
1742 
1743 } /* end mrqpsfexp() */
1744 
1770 void
1771 mrqpsfexp2(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
1772  cxdouble dyda[], cxint na)
1773 {
1774 
1775  const cxchar *fctid = "mrqpsfexp2";
1776 
1777  cxdouble t1,t2,t3,t4,t5,t6,t8,t10,t16;
1778 
1779  cxdouble amplitude = a[LMP_AMPL];
1780  cxdouble center = a[LMP_CENT];
1781  cxdouble background = a[LMP_BACK];
1782  cxdouble width1 = a[LMP_WID1];
1783  cxdouble width2 = a[LMP_WID2];
1784 
1785  /* check for number of parameters */
1786  if (na != 5) {
1787  cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
1788  return;
1789  }
1790 
1791  *y = 0.0;
1792  if (dyda != NULL) {
1793  dyda[LMP_AMPL] = dyda[LMP_CENT] = dyda[LMP_BACK] = dyda[LMP_WID1] =
1794  dyda[LMP_WID2] = 0.0;
1795  }
1796 
1797  t1 = x[0]-center;
1798 
1799  if (t1 > 0.0) {
1800  t2 = t1;
1801  t10 = 1.0;
1802  } else {
1803  t2 = -t1;
1804  t10 = -1.0;
1805  }
1806 
1807  t3 = 1.0/width1;
1808  t4 = t2*t3;
1809  t5 = pow(t4,width2);
1810  t6 = exp(-t5);
1811  t8 = amplitude*t5;
1812  t16 = log(t4);
1813 
1814  *y = amplitude*t6+background;
1815 
1816  /* Check if derivatives expected */
1817  if (dyda == NULL)
1818  return;
1819 
1820  dyda[LMP_AMPL] = t6; /* d(y)/d(amplitude) */
1821  dyda[LMP_CENT] = t8*width2*t10/t2*t6; /* d(y)/d(center) */
1822 
1823  if (isnan(dyda[LMP_CENT]))
1824  dyda[LMP_CENT] = 0.0;
1825 
1826  dyda[LMP_BACK] = 1.0; /* d(y)/d(background) */
1827  dyda[LMP_WID1] = t8*width2*t3*t6; /* d(y)/d(width1) */
1828  dyda[LMP_WID2] = -t8*t16*t6; /* d(y)/d(width2) */
1829 
1830  if (isnan(dyda[LMP_WID2]))
1831  dyda[LMP_WID2] = 0.0;
1832 
1833  if (r != NULL) {
1834  register cxint k;
1835 
1836  k = LMP_AMPL << 1;
1837  if (r[k+1] > 0) {
1838  dyda[LMP_AMPL] *= mrqdydaweight(a[LMP_AMPL],r[k],r[k+1]);
1839  }
1840  k = LMP_CENT << 1;
1841  if (r[k+1] > 0) {
1842  dyda[LMP_CENT] *= mrqdydaweight(a[LMP_CENT],r[k],r[k+1]);
1843  }
1844  k = LMP_WID1 << 1;
1845  if (r[k+1] > 0) {
1846  dyda[LMP_WID1] *= mrqdydaweight(a[LMP_WID1],r[k],r[k+1]);
1847  }
1848  k = LMP_WID2 << 1;
1849  if (r[k+1] > 0) {
1850  dyda[LMP_WID2] *= mrqdydaweight(a[LMP_WID2],r[k],r[k+1]);
1851  }
1852  }
1853 
1854  return;
1855 
1856 } /* end mrqpsfexp2() */
1857 
1878 void
1879 mrqlocywarp(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
1880  cxdouble dyda[], cxint na)
1881 {
1882 
1883  const cxchar *fctid = "mrqlocywarp";
1884 
1885  cxdouble xccd, nx, startx;
1886  register cxint i,ncoef;
1887  cxdouble Tx, Ty, cx, Ky, tt, xx;
1888  cpl_matrix *mBase = NULL, *mX = NULL;
1889  register cxdouble fxx = 0.0, f1xx = 0.0, f2xx = 0.0, z1;
1890  cxdouble *pd_x = NULL, *pd_mbase = NULL;
1891 
1892  /* check for number of parameters */
1893  if (na != 5) {
1894  cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
1895  return;
1896  }
1897 
1898  *y = 0.0;
1899  if (dyda != NULL) {
1900  dyda[LMP_TX] = dyda[LMP_TY] = dyda[LMP_CX] = dyda[LMP_KY] =
1901  dyda[LMP_TT] = 0.0;
1902  }
1903 
1904  xccd = x[LMI_XCCD]; /* pixel abcissa */
1905  nx = x[LMI_NX]; /* number of pixels in dispersion dir. */
1906  startx = x[LMI_STRX]; /* 1st pixel position */
1907  ncoef = (cxint) x[LMI_NCOF]; /* number of chebyshev coef */
1908 
1909  Tx = a[LMP_TX]; /* Translation X */
1910  Ty = a[LMP_TY]; /* Translation Y */
1911  cx = a[LMP_CX]; /* Scaling X: cx = 1/Kx */
1912  Ky = a[LMP_KY]; /* Scaling Y */
1913  tt = a[LMP_TT]; /* Rotation theta: tt = tan(theta) */
1914 
1915  xx = cx*(xccd-Tx);
1916 
1917  mX = cpl_matrix_new(1,1);
1918  pd_x = cpl_matrix_get_data(mX);
1919  pd_x[0] = xx;
1920 
1921  mBase = giraffe_chebyshev_base1d(startx, nx, ncoef, mX);
1922 
1923  pd_mbase = cpl_matrix_get_data(mBase);
1924 
1925  for (i = 0; i < ncoef; i++)
1926  fxx += pd_mbase[i] * x[i+4];
1927 
1928  if (ncoef > 1) {
1929  for (i = 0; i < (ncoef - 1); i++)
1930  f1xx += pd_mbase[i] * (i+1)*x[i+5];
1931  }
1932 
1933  if (ncoef > 2) {
1934  for (i = 0; i < (ncoef - 2); i++)
1935  f2xx += pd_mbase[i] * (i+2)*x[i+6];
1936  }
1937 
1938  if (mX!=NULL) { cpl_matrix_delete(mX); mX = NULL; pd_x = NULL; }
1939  if (mBase!=NULL) { cpl_matrix_delete(mBase); mBase = NULL; pd_mbase = NULL; }
1940 
1941  z1 = 1.0 - tt*tt + tt*f1xx;
1942  *y = Ky*(fxx-tt*xx)/z1 + Ty;
1943 
1944  /* Check if derivatives expected */
1945  if (dyda == NULL)
1946  return;
1947 
1948  dyda[LMP_TX] = /* d(y)/d(Tx) */
1949  (cx*Ky/z1)*((tt-f1xx) + tt*f2xx*(fxx-tt*xx)/z1);
1950 
1951  dyda[LMP_TY] = 1.0; /* d(y)/d(Ty) */
1952 
1953  dyda[LMP_CX] = /* d(y)/d(cx) */
1954  (Ky*(xccd-Tx)/z1)*((f1xx-tt) - tt*f2xx*(fxx-tt*xx)/z1);
1955 
1956  dyda[LMP_KY] = (fxx-tt*xx)/z1; /* d(y)/d(Ky) */
1957  dyda[LMP_TT] = /* d(y)/d(tt) */
1958  (Ky/(z1*z1))*(-xx*(1.+tt*tt)+2.*tt*fxx-fxx*f1xx);
1959 
1960  if (r != NULL) {
1961  register cxint k;
1962 
1963  k = LMP_TX << 1;
1964  if (r[k+1] > 0) {
1965  dyda[LMP_TX] *= mrqdydaweight(a[LMP_TX],r[k],r[k+1]);
1966  }
1967  k = LMP_CX << 1;
1968  if (r[k+1] > 0) {
1969  dyda[LMP_CX] *= mrqdydaweight(a[LMP_CX],r[k],r[k+1]);
1970  }
1971  k = LMP_KY << 1;
1972  if (r[k+1] > 0) {
1973  dyda[LMP_KY] *= mrqdydaweight(a[LMP_KY],r[k],r[k+1]);
1974  }
1975  k = LMP_TT << 1;
1976  if (r[k+1] > 0) {
1977  dyda[LMP_TT] *= mrqdydaweight(a[LMP_TT],r[k],r[k+1]);
1978  }
1979  }
1980 
1981  return;
1982 
1983 } /* end mrqlocywarp() */
1984 
2005 void
2006 mrqxoptmodGS(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
2007  cxdouble dyda[], cxint na)
2008 {
2009 
2010  const cxchar *fctid = "mrqxoptmodGS";
2011 
2012  register cxdouble lambda,xfibre,yfibre,pixsize,nx;
2013  /* Optical model parameters */
2014  register cxdouble fcoll,cfact;
2015  /* Grating parameters */
2016  register cxdouble gtheta,gorder,gspace;
2017 
2018  register cxdouble t1,t10,t109,t11,t110,t114,t12,t14,t17,t18,t2,t21,t23,t24;
2019  register cxdouble t25,t26,t27,t28,t29,t3,t30,t31,t35,t40,t43,t49,t5,t51,t52;
2020  register cxdouble t53,t59,t6,t66,t67,t69,t7,t71,t8,t82,t84,t9,t95,t98;
2021 
2022  /* check for number of parameters */
2023  if (na != 7) {
2024  cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
2025  return;
2026  }
2027 
2028  *y = 0.0;
2029  if (dyda != NULL) {
2030  dyda[LMP_NX] = dyda[LMP_PXSIZ] = dyda[LMP_FCOLL] = dyda[LMP_CFACT] =
2031  dyda[LMP_THETA] = dyda[LMP_ORDER] = dyda[LMP_SPACE] = 0.0;
2032  }
2033 
2034  lambda = x[LMI_WLEN]; /* wavelength [mm] */
2035  xfibre = x[LMI_XFIB]; /* Y fibre [mm] */
2036  yfibre = x[LMI_YFIB]; /* Y fibre [mm] */
2037 
2038  nx = a[LMP_NX]; /* CCD size in X [pixels] */
2039  pixsize = a[LMP_PXSIZ]; /* CCD pixel size [mm] */
2040  fcoll = a[LMP_FCOLL]; /* collimator focal length [mm] */
2041  cfact = a[LMP_CFACT]; /* camera magnification factor */
2042  gtheta = a[LMP_THETA]; /* grating angle [radian] */
2043  gorder = a[LMP_ORDER]; /* grating diffraction order */
2044  gspace = a[LMP_SPACE]; /* grating groove spacing [mm] */
2045 
2046  t1 = cfact*fcoll;
2047  t2 = lambda*gorder;
2048  t3 = 1.0/gspace;
2049  t5 = cos(gtheta);
2050  t6 = xfibre*t5;
2051  t7 = xfibre*xfibre;
2052  t8 = yfibre*yfibre;
2053  t9 = fcoll*fcoll;
2054  t10 = t7+t8+t9;
2055  t11 = sqrt(t10);
2056  t12 = 1.0/t11;
2057  t14 = sin(gtheta);
2058  t17 = -t2*t3+t12*t6+fcoll*t14*t12;
2059  t18 = t17*t5;
2060  t21 = t17*t17;
2061  t23 = sqrt(1.0-t8/t10-t21);
2062  t24 = t14*t23;
2063  t25 = t18+t24;
2064  t26 = t17*t14;
2065  t27 = t5*t23;
2066  t28 = -t26+t27;
2067  t29 = 1.0/t28;
2068  t30 = t25*t29;
2069  t31 = 1.0/pixsize;
2070  t35 = pixsize*pixsize;
2071  t40 = t29*t31;
2072  t43 = 1.0/t11/t10;
2073  t49 = -t6*t43*fcoll+t14*t12-t9*t14*t43;
2074  t51 = 1.0/t23;
2075  t52 = t14*t51;
2076  t53 = t10*t10;
2077  t59 = 2.0*t8/t53*fcoll-2.0*t17*t49;
2078  t66 = t1*t25;
2079  t67 = t28*t28;
2080  t69 = 1.0/t67*t31;
2081  t71 = t5*t51;
2082  t82 = -xfibre*t14*t12+fcoll*t5*t12;
2083  t84 = t17*t82;
2084  t95 = lambda*t3;
2085  t98 = t17*lambda*t3;
2086  t109 = gspace*gspace;
2087  t110 = 1.0/t109;
2088  t114 = t2*t110;
2089 
2090  /* takes care of model direction */
2091  if (nx < 0.0)
2092  *y = t1*t30*t31-0.5*nx;
2093  else
2094  *y = -t1*t30*t31+0.5*nx;
2095 
2096  /* Check if derivatives expected */
2097  if (dyda == NULL)
2098  return;
2099 
2100  /* derivatives for each parameters */
2101  dyda[LMP_NX] = 0.5; /* d(y)/d(nx) */
2102  dyda[LMP_PXSIZ] = -t1*t30/t35; /* d(y)/d(pixsize) */
2103  dyda[LMP_FCOLL] = /* d(y)/d(fcoll) */
2104  cfact*t25*t40+t1*(t49*t5+t52*t59/2.0)*t29*t31-
2105  t66*t69*(-t49*t14+t71*t59/2.0);
2106  dyda[LMP_CFACT] = /* d(y)/d(cfact) */
2107  fcoll*t25*t40;
2108  dyda[LMP_THETA] = /* d(y)/d(gtheta) */
2109  t1*(t82*t5-t26+t27-t52*t84)*t29*t31-
2110  t66*t69*(-t82*t14-t18-t24-t71*t84);
2111  dyda[LMP_ORDER] = /* d(y)/d(gorder) */
2112  t1*(-t95*t5+t52*t98)*t29*t31-t66*t69*(t95*t14+t71*t98);
2113  dyda[LMP_SPACE] = /* d(y)/d(gspace) */
2114  t1*(t2*t110*t5-t52*t17*t114)*t29*t31-
2115  t66*t69*(-t2*t110*t14-t71*t17*t114);
2116 
2117  if (nx > 0.0) {
2118  dyda[LMP_NX] = -dyda[LMP_NX];
2119  dyda[LMP_PXSIZ] = -dyda[LMP_PXSIZ];
2120  dyda[LMP_FCOLL] = -dyda[LMP_FCOLL];
2121  dyda[LMP_CFACT] = -dyda[LMP_CFACT];
2122  dyda[LMP_THETA] = -dyda[LMP_THETA];
2123  dyda[LMP_ORDER] = -dyda[LMP_ORDER];
2124  dyda[LMP_SPACE] = -dyda[LMP_SPACE];
2125  }
2126 
2127  if (r != NULL) {
2128  register cxint k;
2129 
2130  k = LMP_PXSIZ << 1;
2131  if (r[k+1] > 0) {
2132  dyda[LMP_PXSIZ] *= mrqdydaweight(a[LMP_PXSIZ],r[k],r[k+1]);
2133  }
2134  k = LMP_FCOLL << 1;
2135  if (r[k+1] > 0) {
2136  dyda[LMP_FCOLL] *= mrqdydaweight(a[LMP_FCOLL],r[k],r[k+1]);
2137  }
2138  k = LMP_CFACT << 1;
2139  if (r[k+1] > 0) {
2140  dyda[LMP_CFACT] *= mrqdydaweight(a[LMP_CFACT],r[k],r[k+1]);
2141  }
2142  k = LMP_THETA << 1;
2143  if (r[k+1] > 0) {
2144  dyda[LMP_THETA] *= mrqdydaweight(a[LMP_THETA],r[k],r[k+1]);
2145  }
2146  k = LMP_ORDER << 1;
2147  if (r[k+1] > 0) {
2148  dyda[LMP_ORDER] *= mrqdydaweight(a[LMP_ORDER],r[k],r[k+1]);
2149  }
2150  k = LMP_SPACE << 1;
2151  if (r[k+1] > 0) {
2152  dyda[LMP_SPACE] *= mrqdydaweight(a[LMP_SPACE],r[k],r[k+1]);
2153  }
2154  }
2155 
2156 } /* end mrqxoptmodGS() */
2157 
2181 void
2182 mrqtest(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
2183  cxdouble dyda[], cxint na)
2184 {
2185 
2186  const cxchar *fctid = "mrqtest";
2187 
2188  cxdouble a1 = a[0];
2189  cxdouble b1 = a[1];
2190 
2191  /* check for number of parameters */
2192  if (na != 2) {
2193  cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
2194  return;
2195  }
2196 
2197  *y = 0.0;
2198  *y = a1 * x[0] + b1;
2199 
2200  /* Check if derivatives expected */
2201  if (dyda == NULL)
2202  return;
2203 
2204  dyda[0] = 0.0;
2205  dyda[1] = 0.0;
2206 
2207  return;
2208 
2209 } /* end mrqtest() */
2210 

This file is part of the GIRAFFE Pipeline Reference Manual 2.12.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Mon Mar 24 2014 11:43:52 by doxygen 1.8.2 written by Dimitri van Heesch, © 1997-2004