GIRAFFE Pipeline Reference Manual

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

This file is part of the GIRAFFE Pipeline Reference Manual 2.16.10.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Thu Dec 15 2022 21:18:51 by doxygen 1.9.1 written by Dimitri van Heesch, © 1997-2004