GIRAFFE Pipeline Reference Manual

gimath_lm.c

00001 /* $Id: gimath_lm.c,v 1.6 2006/03/10 15:56:10 rpalsa Exp $
00002  *
00003  * This file is part of the GIRAFFE Pipeline
00004  * Copyright (C) 2002-2006 European Southern Observatory
00005  *
00006  * This program is free software; you can redistribute it and/or modify
00007  * it under the terms of the GNU General Public License as published by
00008  * the Free Software Foundation; either version 2 of the License, or
00009  * (at your option) any later version.
00010  *
00011  * This program is distributed in the hope that it will be useful,
00012  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  * GNU General Public License for more details.
00015  *
00016  * You should have received a copy of the GNU General Public License
00017  * along with this program; if not, write to the Free Software
00018  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00019  */
00020 
00021 /*
00022  * $Author: rpalsa $
00023  * $Date: 2006/03/10 15:56:10 $
00024  * $Revision: 1.6 $
00025  * $Name: giraffe-2_8_8 $
00026  */
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #  include <config.h>
00030 #endif
00031 
00032 #include <stdlib.h>
00033 #include <stdio.h>
00034 #include <math.h>
00035 
00036 #include <cxmemory.h>
00037 #include <cxmessages.h>
00038 
00039 #include "gimacros.h"
00040 #include "gimath.h"
00041 #include "gimath_lm.h"
00042 #include "gimessages.h"
00043 
00058 #define SWAP(a,b) {swap=(a);(a)=(b);(b)=swap;}
00059 
00060 
00068 lmrq_model lmrq_models[] = {
00069     {LMRQ_GAUSSUM,   mrqgaussum,   4, 1, "gaussum",   LINE_MODEL},
00070     {LMRQ_XOPTMOD,   mrqxoptmod,   7, 3, "xoptmod",   XOPT_MODEL},
00071     {LMRQ_XOPTMODGS, mrqxoptmodGS, 7, 3, "xoptmodGS", XOPT_MODEL},
00072     {LMRQ_XOPTMOD2,  mrqxoptmod2, 10, 3, "xoptmod2",  XOPT_MODEL},
00073     {LMRQ_PSFCOS,    mrqpsfcos,    5, 1, "psfcos",    LINE_MODEL},
00074     {LMRQ_PSFEXP,    mrqpsfexp,    5, 1, "psfexp",    LINE_MODEL},
00075     {LMRQ_YOPTMOD,   mrqyoptmod,   7, 3, "yoptmod",   YOPT_MODEL},
00076     {LMRQ_YOPTMOD2,  mrqyoptmod2, 10, 3, "yoptmod2",  YOPT_MODEL},
00077     {LMRQ_LOCYWARP,  mrqlocywarp,  5, 4, "locywarp",  LOCY_MODEL},
00078     {LMRQ_PSFEXP2,   mrqpsfexp2,   5, 1, "psfexp2",   LINE_MODEL},
00079     {LMRQ_TEST,      mrqtest,      2, 1, "test",      LINE_MODEL}
00080 };
00081 
00082 cxint nr_lmrq_models = CX_N_ELEMENTS(lmrq_models);
00083 
00084 
00101 static void
00102 covariance_sort(cpl_matrix *covar, cxint ma, cxint ia[], cxint mfit)
00103 {
00104 
00105     register cxint i, j, k;
00106     register cxdouble swap;
00107 
00108     cxdouble *pd_covar = NULL;
00109     cxint     nr_covar;
00110 
00111     pd_covar = cpl_matrix_get_data(covar);
00112     nr_covar = cpl_matrix_get_nrow(covar);
00113 
00114     for (i = mfit; i < ma; i++)
00115         for (j = 0; j <= i; j++)
00116             pd_covar[i * nr_covar + j] = pd_covar[j * nr_covar + i] = 0.0;
00117 
00118     k = mfit - 1;
00119     for (j = (ma-1); j >= 0; j--) {
00120         if (ia[j]) {
00121             for (i = 0; i < ma; i++)
00122                 SWAP(pd_covar[i * nr_covar + k],pd_covar[i * nr_covar + j])
00123 
00124             for (i = 0;i < ma; i++)
00125                 SWAP(pd_covar[k * nr_covar + i],pd_covar[j * nr_covar + i])
00126 
00127             k--;
00128         }
00129     }
00130 
00131 } /* end covariance_sort() */
00132 
00133 #undef SWAP
00134 
00161 static double
00162 mrqdydaweight(cxdouble x, cxdouble x0, cxdouble dx)
00163 {
00164     register cxdouble w;
00165 
00166     w = exp(-pow(fabs(x-x0),DW_DEGREE)/pow(dx,DW_DEGREE/DW_LOG001));
00167 
00168     if (isnan(w))
00169         w = 1;
00170 
00171     return w;
00172 }
00173 
00200 cxint
00201 mrqnlfit(
00202     cpl_matrix  *x,
00203     cpl_matrix  *y,
00204     cpl_matrix  *sig,
00205     cxint        ndata,
00206     cpl_matrix  *a,
00207     cxdouble     r[],
00208     cxint        ia[],
00209     cxint        ma,
00210     cpl_matrix  *alpha,
00211     cxdouble    *chisq,
00212     lmrq_params  fit_params,
00213     fitted_func  funcs
00214 ) {
00215 
00216     cxint       itst,
00217                 n,
00218                 res;
00219 
00220     cxdouble    alamda,
00221                 ochisq;
00222 
00223     cpl_matrix *beta = NULL;
00224 
00225     /*************************************************************************
00226                                       PROCESSING
00227     *************************************************************************/
00228 
00229     beta = cpl_matrix_new(ma,ma);
00230 
00231     alamda = -1.0;
00232 
00233     res = mymrqmin(x, y, sig, ndata, a, r, ia, ma, alpha, beta, chisq,
00234                    funcs, &alamda);
00235 
00236     if (res != 0) {
00237         cpl_matrix_delete(beta); beta = NULL;
00238         return res;
00239     }
00240 
00241     itst=0;
00242 
00243     for (n = 1; n <= fit_params.imax; n++) {
00244 
00245         ochisq = *chisq;
00246 
00247         res = mymrqmin(x, y, sig, ndata, a, r, ia, ma, alpha, beta, chisq,
00248                        funcs, &alamda);
00249 
00250         if (res!=0) {
00251             cpl_matrix_delete(beta); beta = NULL;
00252             return res;
00253         }
00254 
00255         if (*chisq > ochisq)
00256             itst=0;
00257         else if (fabs(ochisq-*chisq) < fit_params.dchsq)
00258             itst++;
00259 
00260         if (itst > fit_params.tmax)
00261             break;
00262     }
00263 
00264     /* get covariance matrix */
00265     alamda=0.0;
00266 
00267     res = mymrqmin(x, y, sig, ndata, a, r, ia, ma, alpha, beta, chisq,
00268                    funcs, &alamda);
00269 
00270     if (res != 0) {
00271         cpl_matrix_delete(beta); beta = NULL;
00272         return res;
00273     }
00274 
00275     cpl_matrix_delete(beta); beta = NULL;
00276 
00277     return n;
00278 
00279 } /* end mrqnlfit() */
00280 
00337 cxint
00338 mymrqmin(
00339     cpl_matrix  *x,
00340     cpl_matrix  *y,
00341     cpl_matrix  *sig,
00342     cxint        ndata,
00343     cpl_matrix  *a,
00344     cxdouble     r[],
00345     cxint        ia[],
00346     cxint        ma,
00347     cpl_matrix  *covar,
00348     cpl_matrix  *alpha,
00349     cxdouble    *chisq,
00350     fitted_func  funcs,
00351     cxdouble    *alamda
00352 ) {
00353 
00354     register cxint    gj, j, k, l, m;
00355 
00356     static cxdouble   *pd_a, *pd_covar, *pd_alpha;
00357     static cxint       nr_covar, nr_alpha, nr_moneda, mfit;
00358 
00359     static cpl_matrix *matry, *mbeta, *mda, *moneda;
00360     static cxdouble   *atry,  *beta,  *da,  *oneda, ochisq;
00361 
00362     /*************************************************************************
00363                                       PROCESSING
00364     *************************************************************************/
00365 
00366     pd_a     = cpl_matrix_get_data(a);
00367     pd_covar = cpl_matrix_get_data(covar);
00368     pd_alpha = cpl_matrix_get_data(alpha);
00369     nr_covar = cpl_matrix_get_nrow(covar);
00370     nr_alpha = cpl_matrix_get_nrow(alpha);
00371 
00372     if (*alamda<0.0) {
00373 
00374         matry = cpl_matrix_new(ma,1);
00375         atry  = cpl_matrix_get_data(matry);
00376 
00377         mbeta = cpl_matrix_new(ma,1);
00378         beta  = cpl_matrix_get_data(mbeta);
00379 
00380         mda = cpl_matrix_new(ma,1);
00381         da  = cpl_matrix_get_data(mda);
00382 
00383         for (mfit = 0, j = 0; j < ma; j++)
00384             if (ia[j])
00385                 mfit++;
00386 
00387         moneda    = cpl_matrix_new(1,mfit);
00388         oneda     = cpl_matrix_get_data(moneda);
00389 
00390         *alamda = 0.001;
00391 
00392         gj = mymrqcof(x, y, sig, ndata, a, r, ia, ma, alpha, mbeta,
00393                       chisq, funcs);
00394 
00395         if (gj != 0) {
00396             cpl_matrix_delete(moneda); moneda = NULL; oneda = NULL;
00397             cpl_matrix_delete(mda);    mda    = NULL; da    = NULL;
00398             cpl_matrix_delete(mbeta);  mbeta  = NULL; beta  = NULL;
00399             cpl_matrix_delete(matry);  matry  = NULL; atry  = NULL;
00400             return gj;
00401         }
00402 
00403         ochisq = (*chisq);
00404 
00405         for (j = 0; j < ma; j++)
00406             atry[j] = pd_a[j];
00407 
00408     }
00409 
00410     nr_moneda = cpl_matrix_get_nrow(moneda);
00411 
00412     for (j = -1, l = 0; l < ma; l++) {
00413         if (ia[l]) {
00414             for (j++, k = -1, m = 0; m < ma; m++) {
00415                 if (ia[m]) {
00416                     k++;
00417                     pd_covar[j * nr_covar + k] = pd_alpha[j * nr_alpha + k];
00418                 }
00419             }
00420 
00421             pd_covar[j * nr_covar + j] =
00422                 pd_alpha[j * nr_alpha + j] * (1.0 + (*alamda));
00423 
00424             oneda[j * nr_moneda + 0] = beta[j];
00425         }
00426     }
00427 
00428     gj = giraffe_gauss_jordan(covar, mfit, moneda, 1);
00429 
00430     if (gj != 0) {
00431         cpl_matrix_delete(moneda); moneda = NULL; oneda = NULL;
00432         cpl_matrix_delete(mda);    mda    = NULL; da    = NULL;
00433         cpl_matrix_delete(mbeta);  mbeta  = NULL; beta  = NULL;
00434         cpl_matrix_delete(matry);  matry  = NULL; atry  = NULL;
00435         return gj;
00436     }
00437 
00438     for (j = 0; j < mfit; j++)
00439         da[j] = oneda[j * nr_moneda + 0];
00440 
00441     if (*alamda == 0.0) {
00442         covariance_sort(covar, ma, ia, mfit);
00443         cpl_matrix_delete(moneda); moneda = NULL; oneda = NULL;
00444         cpl_matrix_delete(mda);    mda    = NULL; da    = NULL;
00445         cpl_matrix_delete(mbeta);  mbeta  = NULL; beta  = NULL;
00446         cpl_matrix_delete(matry);  matry  = NULL; atry  = NULL;
00447         return 0;
00448     }
00449 
00450     for (j = -1, l = 0; l < ma; l++)
00451         if (ia[l])
00452             atry[l] = pd_a[l] + da[++j];
00453 
00454     gj = mymrqcof(x, y, sig, ndata, matry, r, ia, ma, covar, mda,
00455                   chisq, funcs);
00456 
00457     if (gj != 0) {
00458         cpl_matrix_delete(moneda); moneda = NULL; oneda = NULL;
00459         cpl_matrix_delete(mda);    mda    = NULL; da    = NULL;
00460         cpl_matrix_delete(mbeta);  mbeta  = NULL; beta  = NULL;
00461         cpl_matrix_delete(matry);  matry  = NULL; atry  = NULL;
00462         return gj;
00463     }
00464 
00465     if (*chisq < ochisq) {
00466 
00467         *alamda *= 0.1;
00468         ochisq = *chisq;
00469 
00470         for (j = -1, l = 0; l < ma; l++) {
00471             if (ia[l]) {
00472                 for (j++, k = -1, m = 0; m < ma; m++) {
00473                     if (ia[m]) {
00474                         k++;
00475                         pd_alpha[j * nr_alpha + k] =
00476                             pd_covar[j * nr_covar + k];
00477                     }
00478                 }
00479 
00480                 beta[j] = da[j];
00481                 pd_a[l] = atry[l];
00482             }
00483         }
00484 
00485     } else {
00486         *alamda *= 10.0;
00487         *chisq = ochisq;
00488     }
00489 
00490     return 0;
00491 
00492 } /* end mymrqmin() */
00493 
00521 cxint
00522 mymrqcof(
00523     cpl_matrix  *x,
00524     cpl_matrix  *y,
00525     cpl_matrix  *sig,
00526     cxint        ndata,
00527     cpl_matrix  *a,
00528     cxdouble     r[],
00529     cxint        ia[],
00530     cxint        ma,
00531     cpl_matrix  *alpha,
00532     cpl_matrix  *beta,
00533     cxdouble    *chisq,
00534     fitted_func  funcs
00535 ) {
00536 
00537     register cxint i, j, k, l, m, mfit = 0;
00538 
00539     cxdouble  ymod, wt, sig2i, dy, *dyda;
00540 
00541     cxdouble *pd_x     = NULL,
00542              *pd_y     = NULL,
00543              *pd_sig   = NULL,
00544              *pd_a     = NULL,
00545              *pd_alpha = NULL,
00546              *pd_beta  = NULL;
00547 
00548     cxint     nr_alpha, nc_x;
00549 
00550     /************************************************************************
00551                                      PROCESSING
00552     ************************************************************************/
00553 
00554     pd_x     = cpl_matrix_get_data(x);
00555     nc_x     = cpl_matrix_get_ncol(x);
00556     pd_y     = cpl_matrix_get_data(y);
00557     pd_sig   = cpl_matrix_get_data(sig);
00558     pd_a     = cpl_matrix_get_data(a);
00559     pd_alpha = cpl_matrix_get_data(alpha);
00560     nr_alpha = cpl_matrix_get_nrow(alpha);
00561     pd_beta  = cpl_matrix_get_data(beta);
00562 
00563     for (j = 0; j < ma; j++) {
00564         if (ia[j])
00565             mfit++;
00566     }
00567 
00568     for (j = 0; j < mfit; j++) {
00569         for (k = 0; k <= j; k++)
00570             pd_alpha[j * nr_alpha + k] = 0.0;
00571 
00572         pd_beta[j] = 0.0;
00573     }
00574 
00575     *chisq = 0.0;
00576 
00577     dyda = (cxdouble *) cx_calloc(ma, sizeof(cxdouble));
00578 
00579     for (i = 0; i < ndata; i++) {
00580 
00581         (*funcs)(&(pd_x[i*nc_x]), pd_a, r, &ymod, dyda, ma);
00582 
00583         if (pd_sig[i]==0.0) {
00584             continue;
00585         }
00586 
00587         sig2i = 1.0 / (pd_sig[i] * pd_sig[i]);
00588 
00589         dy = pd_y[i] - ymod;
00590 
00591         for (j = -1, l = 0; l < ma; l++) {
00592 
00593            if (ia[l]) {
00594                 wt = dyda[l] * sig2i;
00595                 for (j++, k = -1, m = 0; m <= l; m++) {
00596                     if (ia[m]) {
00597                         ++k;
00598                         pd_alpha[j * nr_alpha + k] += (wt * dyda[m]);
00599                     }
00600                 }
00601 
00602                 pd_beta[j] += (dy * wt);
00603 
00604             }
00605         }
00606 
00607         *chisq += (dy * dy * sig2i);
00608 
00609     }
00610 
00611     for (j = 1; j < mfit; j++)
00612         for (k = 0; k < j; k++)
00613             pd_alpha[k * nr_alpha + j] = pd_alpha[j * nr_alpha + k];
00614 
00615 
00616     cx_free(dyda);
00617 
00618     return 0;
00619 
00620 } /* end mymrqcof() */
00621 
00638 cxdouble
00639 r_squared(cxdouble resSS, cpl_matrix *y, cxint n)
00640 {
00641     register cxint i;
00642     register cxdouble Sy, Syy, SS;
00643     cxdouble res, *pd_y = NULL;
00644 
00645     pd_y = cpl_matrix_get_data(y);
00646 
00647     if (n < 1)
00648         return 0.0;
00649 
00650     for (i=0, Sy=0.0, Syy=0.0; i<n; i++) {
00651         Sy  += pd_y[i];
00652         Syy += pd_y[i]*pd_y[i];
00653     }
00654 
00655     SS = Syy - Sy*Sy/n;
00656     res = resSS/SS;
00657 
00658     if (isnan(res))
00659         return 0.0;
00660 
00661     if (res > 0.0)
00662         res = sqrt(res);
00663 
00664     return res;
00665 
00666 } /* end r_squared() */
00667 
00699 void
00700 mrqgaussum(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
00701            cxdouble dyda[], cxint na)
00702 {
00703 
00704     register cxint i, j;
00705     register cxdouble fac,ex,amplitude,center,backg,width,xred;
00706 
00707     *y = 0.0;
00708 
00709     for (j = 0, i = 0; i < na; i += 4, j += 4) {
00710         amplitude = a[i];
00711         center    = a[i + 1];
00712         backg     = a[i + 2];
00713         width     = a[i + 3];
00714         xred = (x[0] - center) / width;
00715         ex   = exp(-xred * xred / 2.);
00716         fac  = amplitude * xred * ex;
00717         *y  += (amplitude * ex + backg);
00718 
00719         /* Check if derivatives expected */
00720         if (dyda == NULL) continue;
00721 
00722         /* derivatives for each parameters */
00723         dyda[j]     = ex;                      /* d(y)/d(amplitude) */
00724         dyda[j + 1] = fac / width;             /* d(y)/d(center)    */
00725         dyda[j + 2] = 1.;                      /* d(y)/d(backg)     */
00726         dyda[j + 3] = (fac * xred) / width;    /* d(y)/d(width)     */
00727     }
00728 
00729 } /* end mrqgaussum() */
00730 
00782 void
00783 mrqxoptmod(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
00784            cxdouble dyda[], cxint na)
00785 {
00786 
00787     const cxchar *fctid = "mrqxoptmod";
00788 
00789     register cxdouble xccd, d, X;
00790     register cxdouble lambda,xfibre,yfibre,pixsize,nx;
00791     /* Optical model parameters */
00792     register cxdouble fcoll,cfact;
00793     /* Grating parameters */
00794     register cxdouble gtheta,gorder,gspace;
00795     register cxdouble yfibre2,tmp,tmp2,d2,X2,gspace2,sqtmp,costheta,sintheta;
00796 
00797     /* check for number of parameters */
00798     if (na != 7) {
00799         cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
00800         return;
00801     }
00802 
00803     *y = 0.0;
00804     if (dyda != NULL) {
00805         dyda[0] = dyda[1] = dyda[2] = dyda[3] =
00806                       dyda[4] = dyda[5] = dyda[6] = 0.0;
00807     }
00808 
00809     lambda  = x[LMI_WLEN];    /* wavelength [mm]              */
00810     xfibre  = x[LMI_XFIB];    /* Y fibre [mm]                 */
00811     yfibre  = x[LMI_YFIB];    /* Y fibre [mm]                 */
00812 
00813     nx      = a[LMP_NX];      /* CCD  size in X [pixels]      */
00814     pixsize = a[LMP_PXSIZ];   /* CCD pixel size [mm]          */
00815     fcoll   = a[LMP_FCOLL];   /* collimator focal length [mm] */
00816     cfact   = a[LMP_CFACT];   /* camera magnification factor  */
00817     gtheta  = a[LMP_THETA];   /* grating angle [radian]       */
00818     gorder  = a[LMP_ORDER];   /* grating diffraction order    */
00819     gspace  = a[LMP_SPACE];   /* grating groove spacing [mm]  */
00820 
00821     yfibre2  = yfibre * yfibre;
00822     gspace2  = gspace * gspace;
00823     costheta = cos(gtheta);
00824     sintheta = sin(gtheta);
00825     d2 = xfibre * xfibre + yfibre2 + (fcoll * fcoll);
00826     d  = sqrt(d2);
00827     X  = (-lambda*gorder/gspace) + (xfibre*costheta/d) + (fcoll*sintheta/d);
00828     X2 = X * X;
00829     sqtmp = sqrt(1.0 - yfibre2/d2 - X2);
00830     tmp   = -sintheta*X + costheta*sqtmp;
00831     tmp2  = tmp * tmp;
00832     xccd  = (cfact * fcoll * (X*costheta + sintheta*sqtmp))/tmp;
00833 
00834     /* takes care of model direction */
00835     if (nx < 0.0)
00836         *y = (xccd / pixsize - 0.5*nx);
00837     else
00838         *y = (-xccd / pixsize + 0.5*nx);
00839 
00840     /* Check if derivatives expected */
00841     if (dyda == NULL)
00842         return;
00843 
00844     /* derivatives for each parameters */
00845     dyda[LMP_NX]    = 0.5;                   /* d(y)/d(nx)      */
00846     dyda[LMP_PXSIZ] = 0.0;                   /* d(y)/d(pixsize) */
00847 
00848     dyda[LMP_FCOLL] = cfact*(costheta*X+sintheta*sqtmp)/tmp +
00849               cfact*fcoll*(costheta*(-X*fcoll/d2+sintheta/d -
00850                gorder*lambda*fcoll/(d2*gspace)) +
00851                0.5*sintheta*(-2.0*X*(-X*fcoll/d2+sintheta/d -
00852                gorder*lambda*fcoll/(d2*gspace))+2.0*yfibre2*fcoll/(d2*d2))/sqtmp)/tmp -
00853               cfact*fcoll*(costheta*X+sintheta*sqtmp)*(-sintheta*(-X*fcoll/d2 +
00854                sintheta/d-gorder*lambda*fcoll/(d2*gspace)) +
00855                0.5*costheta*(-2.0*X*(-X*fcoll/d2+sintheta/d -
00856                gorder*lambda*fcoll/(d2*gspace))+2.0*yfibre2*fcoll/(d2*d2))/sqtmp)/tmp2;
00857     dyda[LMP_FCOLL] /= pixsize;                 /* d(y)/d(fcoll) */
00858 
00859     dyda[LMP_CFACT] = (xccd/cfact)/pixsize;     /* d(y)/d(cfact) */
00860 
00861     dyda[LMP_THETA] = cfact*fcoll*((-xfibre*sintheta/d+fcoll*costheta/d)*costheta -
00862                sintheta*X-sintheta*X*(-xfibre*sintheta/d+fcoll*costheta/d)/sqtmp +
00863                costheta*sqtmp)/tmp -
00864               cfact*fcoll*(costheta*X+sintheta*sqtmp)*(-(-xfibre*sintheta/d +
00865                fcoll*costheta/d)*sintheta-costheta*X -
00866                costheta*X*(-xfibre*sintheta/d+fcoll*costheta/d)/sqtmp -
00867                sintheta*sqtmp)/tmp2;
00868     dyda[LMP_THETA] /= pixsize;             /* d(y)/d(gtheta) */
00869 
00870     dyda[LMP_ORDER] = 0.0;                   /* d(y)/d(gorder) */
00871     dyda[LMP_SPACE] = cfact*fcoll*(lambda*gorder*costheta/gspace2-sintheta*X*lambda*gorder/(sqtmp*gspace2))/tmp -
00872               cfact*fcoll*(X*costheta+sintheta*sqtmp) *
00873               (-lambda*gorder*sintheta/gspace2-costheta*X*lambda*gorder/(sqtmp*gspace2))/tmp2;
00874     dyda[LMP_SPACE] /= pixsize;             /* d(y)/d(gspace) */
00875 
00876     if (nx > 0.0) {
00877         dyda[LMP_NX]    = -dyda[LMP_NX];
00878         dyda[LMP_PXSIZ] = -dyda[LMP_PXSIZ];
00879         dyda[LMP_FCOLL] = -dyda[LMP_FCOLL];
00880         dyda[LMP_CFACT] = -dyda[LMP_CFACT];
00881         dyda[LMP_THETA] = -dyda[LMP_THETA];
00882         dyda[LMP_ORDER] = -dyda[LMP_ORDER];
00883         dyda[LMP_SPACE] = -dyda[LMP_SPACE];
00884     }
00885 
00886     if (r != NULL) {
00887         register cxint k;
00888 
00889         k = LMP_FCOLL << 1;
00890         if (r[k+1] > 0) {
00891             dyda[LMP_FCOLL] *= mrqdydaweight(a[LMP_FCOLL],r[k],r[k+1]);
00892         }
00893         k = LMP_CFACT << 1;
00894         if (r[k+1] > 0) {
00895             dyda[LMP_CFACT] *= mrqdydaweight(a[LMP_CFACT],r[k],r[k+1]);
00896         }
00897         k = LMP_THETA << 1;
00898         if (r[k+1] > 0) {
00899             dyda[LMP_THETA] *=  mrqdydaweight(a[LMP_THETA],r[k],r[k+1]);
00900         }
00901         k = LMP_SPACE << 1;
00902         if (r[k+1] > 0) {
00903             dyda[LMP_SPACE] *=  mrqdydaweight(a[LMP_SPACE],r[k],r[k+1]);
00904         }
00905     }
00906 
00907 } /* end mrqxoptmod() */
00908 
00970 void
00971 mrqxoptmod2(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
00972             cxdouble dyda[], cxint na)
00973 {
00974 
00975     const cxchar *fctid = "mrqxoptmod2";
00976 
00977     register cxdouble lambda,xfibre,yfibre,pixsize,nx;
00978     /* Optical model parameters */
00979     register cxdouble fcoll,cfact;
00980     /* Grating parameters */
00981     register cxdouble gtheta,gorder,gspace;
00982     /* Slit position parameters */
00983     cxdouble  slitdx,slitdy,slitphi;
00984 
00985     register cxdouble t1,t10,t104,t107,t11,t113,t119,t12,t120,t121,t124,t136,
00986                     t137,t138,t14,t143,t148,t16,t161,t162,t166,t168,t17,t173,
00987                     t18,t19,t191,t195,t196,t2,t20,t201,t21,t210,t23,t24,t26,
00988                     t27,t28,t3,t30,t32,t33,t34,t35,t36,t37,t38,t39,t4,t40,t44,
00989                     t49,t52,t58,t60,t61,t62,t64,t68,t75,t76,t78,t80,t9,t91,t93;
00990 
00991     /* check for number of parameters */
00992     if (na != 10) {
00993         cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
00994         return;
00995     }
00996 
00997     *y = 0.0;
00998     if (dyda != NULL) {
00999         dyda[0] = dyda[1] = dyda[2] = dyda[3] =
01000                       dyda[4] = dyda[5] = dyda[6] =
01001                       dyda[7] = dyda[8] = dyda[9] = 0.0;
01002     }
01003 
01004     lambda  = x[LMI_WLEN];    /* wavelength [mm]              */
01005     xfibre  = x[LMI_XFIB];    /* Y fibre [mm]                 */
01006     yfibre  = x[LMI_YFIB];    /* Y fibre [mm]                 */
01007 
01008     nx      = a[LMP_NX];      /* CCD  size in X [pixels]      */
01009     pixsize = a[LMP_PXSIZ];   /* CCD pixel size [mm]          */
01010     fcoll   = a[LMP_FCOLL];   /* collimator focal length [mm] */
01011     cfact   = a[LMP_CFACT];   /* camera magnification factor  */
01012     gtheta  = a[LMP_THETA];   /* grating angle [radian]       */
01013     gorder  = a[LMP_ORDER];   /* grating diffraction order    */
01014     gspace  = a[LMP_SPACE];   /* grating groove spacing [mm]  */
01015     slitdx  = a[LMP_SOFFX];    /* slit position x offset [mm]  */
01016     slitdy  = a[LMP_SOFFY];    /* slit position y offset [mm]  */
01017     slitphi = a[LMP_SPHI];     /* slit position angle [radian] */
01018 
01019     t1 = cfact*fcoll;
01020     t2 = cos(gtheta);
01021     t3 = lambda*gorder;
01022     t4 = 1.0/gspace;
01023     t9 = xfibre*(1.0+slitphi*yfibre)+slitdx;
01024     t10 = t9*t2;
01025     t11 = t9*t9;
01026     t12 = slitphi*slitphi;
01027     t14 = sqrt(1.0-t12);
01028     t16 = yfibre*t14+slitdy;
01029     t17 = t16*t16;
01030     t18 = fcoll*fcoll;
01031     t19 = t11+t17+t18;
01032     t20 = sqrt(t19);
01033     t21 = 1.0/t20;
01034     t23 = sin(gtheta);
01035     t24 = fcoll*t23;
01036     t26 = -t3*t4+t10*t21+t24*t21;
01037     t27 = t2*t26;
01038     t28 = 1.0/t19;
01039     t30 = t26*t26;
01040     t32 = sqrt(1.0-t17*t28-t30);
01041     t33 = t23*t32;
01042     t34 = t27+t33;
01043     t35 = t23*t26;
01044     t36 = t2*t32;
01045     t37 = -t35+t36;
01046     t38 = 1.0/t37;
01047     t39 = t34*t38;
01048     t40 = 1.0/pixsize;
01049     t44 = pixsize*pixsize;
01050     t49 = t38*t40;
01051     t52 = 1.0/t20/t19;
01052     t58 = -t10*t52*fcoll+t23*t21-t18*t23*t52;
01053     t60 = 1.0/t32;
01054     t61 = t23*t60;
01055     t62 = t19*t19;
01056     t64 = t17/t62;
01057     t68 = 2.0*t64*fcoll-2.0*t26*t58;
01058     t75 = t1*t34;
01059     t76 = t37*t37;
01060     t78 = 1.0/t76*t40;
01061     t80 = t2*t60;
01062     t91 = -t9*t23*t21+fcoll*t2*t21;
01063     t93 = t26*t91;
01064     t104 = t2*lambda;
01065     t107 = t26*lambda*t4;
01066     t113 = t23*lambda;
01067     t119 = gspace*gspace;
01068     t120 = 1.0/t119;
01069     t121 = gorder*t120;
01070     t124 = t3*t120;
01071     t136 = t2*t21;
01072     t137 = 2.0*t9;
01073     t138 = t52*t137;
01074     t143 = t136-t10*t138/2.0-t24*t138/2.0;
01075     t148 = t64*t137-2.0*t26*t143;
01076     t161 = 2.0*t16;
01077     t162 = t52*t161;
01078     t166 = -t10*t162/2.0-t24*t162/2.0;
01079     t168 = t16*t28;
01080     t173 = -2.0*t168+t64*t161-2.0*t26*t166;
01081     t191 = 1.0/t14;
01082     t195 = 2.0*t9*xfibre*yfibre-2.0*t16*yfibre*t191*slitphi;
01083     t196 = t52*t195;
01084     t201 = xfibre*yfibre*t136-t10*t196/2.0-t24*t196/2.0;
01085     t210 = 2.0*t168*yfibre*t191*slitphi+t64*t195-2.0*t26*t201;
01086 
01087     /* takes care of model direction */
01088     if (nx < 0.0)
01089         *y = t1*t39*t40-0.5*nx;
01090     else
01091         *y = -t1*t39*t40+0.5*nx;
01092 
01093     /* Check if derivatives expected */
01094     if (dyda == NULL)
01095         return;
01096 
01097     /* derivatives for each parameters */
01098     dyda[LMP_NX]    = 0.5;                   /* d(y)/d(nx)      */
01099     dyda[LMP_PXSIZ] = -t1*t39/t44;           /* d(y)/d(pixsize) */
01100     dyda[LMP_FCOLL] =                        /* d(y)/d(fcoll) */
01101                     cfact*t34*t49+t1*(t2*t58+t61*t68/2.0)*t38*t40 -
01102                     t75*t78*(-t23*t58+t80*t68/2.0);
01103     dyda[LMP_CFACT] =                        /* d(y)/d(cfact) */
01104                     fcoll*t34*t49;
01105     dyda[LMP_THETA] =                        /* d(y)/d(gtheta) */
01106                     t1*(-t35+t2*t91+t36-t61*t93)*t38*t40 -
01107                     t75*t78*(-t27-t23*t91-t33-t80*t93);
01108     dyda[LMP_ORDER] =                        /* d(y)/d(gorder) */
01109                     t1*(-t104*t4+t61*t107)*t38*t40-t75*t78*(t113*t4+t80*t107);
01110     dyda[LMP_SPACE] =                        /* d(y)/d(gspace) */
01111                     t1*(t104*t121-t61*t26*t124)*t38*t40 -
01112                     t75*t78*(-t113*t121-t80*t26*t124);
01113     dyda[LMP_SOFFX] =                        /* d(y)/d(slitdx) */
01114                     t1*(t2*t143+t61*t148/2.0)*t38*t40 -
01115                     t75*t78*(-t23*t143+t80*t148/2.0);
01116     dyda[LMP_SOFFY] =                        /* d(y)/d(slitdy) */
01117                     t1*(t2*t166+t61*t173/2.0)*t38*t40 -
01118                     t75*t78*(-t23*t166+t80*t173/2.0);
01119     dyda[LMP_SPHI]  =                        /* d(y)/d(slitphi) */
01120                     t1*(t2*t201+t61*t210/2.0)*t38*t40 -
01121                     t75*t78*(-t23*t201+t80*t210/2.0);
01122 
01123     if (nx > 0.0) {
01124         dyda[LMP_NX]    = -dyda[LMP_NX];
01125         dyda[LMP_PXSIZ] = -dyda[LMP_PXSIZ];
01126         dyda[LMP_FCOLL] = -dyda[LMP_FCOLL];
01127         dyda[LMP_CFACT] = -dyda[LMP_CFACT];
01128         dyda[LMP_THETA] = -dyda[LMP_THETA];
01129         dyda[LMP_ORDER] = -dyda[LMP_ORDER];
01130         dyda[LMP_SPACE] = -dyda[LMP_SPACE];
01131         dyda[LMP_SOFFX] = -dyda[LMP_SOFFX];
01132         dyda[LMP_SOFFY] = -dyda[LMP_SOFFY];
01133         dyda[LMP_SPHI]  = -dyda[LMP_SPHI];
01134     }
01135 
01136     if (r != NULL) {
01137         register cxint k;
01138 
01139         k = LMP_PXSIZ << 1;
01140         if (r[k+1] > 0) {
01141             dyda[LMP_PXSIZ] *= mrqdydaweight(a[LMP_PXSIZ],r[k],r[k+1]);
01142         }
01143         k = LMP_FCOLL << 1;
01144         if (r[k+1] > 0) {
01145             dyda[LMP_FCOLL] *= mrqdydaweight(a[LMP_FCOLL],r[k],r[k+1]);
01146         }
01147         k = LMP_CFACT << 1;
01148         if (r[k+1] > 0) {
01149             dyda[LMP_CFACT] *= mrqdydaweight(a[LMP_CFACT],r[k],r[k+1]);
01150         }
01151         k = LMP_THETA << 1;
01152         if (r[k+1] > 0) {
01153             dyda[LMP_THETA] *=  mrqdydaweight(a[LMP_THETA],r[k],r[k+1]);
01154         }
01155         k = LMP_ORDER << 1;
01156         if (r[k+1] > 0) {
01157             dyda[LMP_ORDER] *=  mrqdydaweight(a[LMP_ORDER],r[k],r[k+1]);
01158         }
01159         k = LMP_SPACE << 1;
01160         if (r[k+1] > 0) {
01161             dyda[LMP_SPACE] *=  mrqdydaweight(a[LMP_SPACE],r[k],r[k+1]);
01162         }
01163         k = LMP_SOFFX << 1;
01164         if (r[k+1] > 0) {
01165             dyda[LMP_SOFFX] *=  mrqdydaweight(a[LMP_SOFFX],r[k],r[k+1]);
01166         }
01167         k = LMP_SOFFY << 1;
01168         if (r[k+1] > 0) {
01169             dyda[LMP_SOFFY] *=  mrqdydaweight(a[LMP_SOFFY],r[k],r[k+1]);
01170         }
01171         k = LMP_SPHI << 1;
01172         if (r[k+1] > 0) {
01173             dyda[LMP_SPHI] *=  mrqdydaweight(a[LMP_SPHI],r[k],r[k+1]);
01174         }
01175     }
01176 
01177 } /* end mrqxoptmod2() */
01178 
01243 void
01244 mrqyoptmod(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
01245            cxdouble dyda[], cxint na)
01246 {
01247 
01248     const cxchar *fctid = "mrqyoptmod";
01249 
01250     register cxdouble lambda,xfibre,yfibre,pixsize,ny;
01251     /* Optical model parameters */
01252     register cxdouble fcoll,cfact;
01253     /* Grating parameters */
01254     register cxdouble gtheta,gorder,gspace;
01255 
01256     cxdouble t10,t12,t13,t15,t18,t2,t22,t24,t26,t27,t28,t29,t3,t30,t33,
01257            t4,t41,t45,t47,t5,t53,t56,t57,t6,t7,t76,t8,t9,t93,t94;
01258 
01259     /* check for number of parameters */
01260     if (na != 7) {
01261         cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
01262         return;
01263     }
01264 
01265     *y = 0.0;
01266     if (dyda != NULL) {
01267         dyda[LMP_NX] = dyda[LMP_PXSIZ] = dyda[LMP_FCOLL] = dyda[LMP_CFACT] =
01268                       dyda[LMP_THETA] = dyda[LMP_ORDER] = dyda[LMP_SPACE] = 0.0;
01269     }
01270 
01271     lambda  = x[LMI_WLEN];
01272     xfibre  = x[LMI_XFIB];
01273     yfibre  = x[LMI_YFIB];
01274 
01275     ny      = a[LMP_NY];
01276     pixsize = a[LMP_PXSIZ];
01277     fcoll   = a[LMP_FCOLL];
01278     cfact   = a[LMP_CFACT];
01279     gtheta  = a[LMP_THETA];
01280     gorder  = a[LMP_ORDER];
01281     gspace  = a[LMP_SPACE];
01282 
01283     t2 = cfact*fcoll*yfibre;
01284     t3 = xfibre*xfibre;
01285     t4 = yfibre*yfibre;
01286     t5 = fcoll*fcoll;
01287     t6 = t3+t4+t5;
01288     t7 = sqrt(t6);
01289     t8 = 1.0/t7;
01290     t9 = lambda*gorder;
01291     t10 = 1.0/gspace;
01292     t12 = cos(gtheta);
01293     t13 = xfibre*t12;
01294     t15 = sin(gtheta);
01295     t18 = -t9*t10+t13*t8+fcoll*t15*t8;
01296     t22 = t18*t18;
01297     t24 = sqrt(1.0-t4/t6-t22);
01298     t26 = -t18*t15+t12*t24;
01299     t27 = 1.0/t26;
01300     t28 = t8*t27;
01301     t29 = 1.0/pixsize;
01302     t30 = t28*t29;
01303     t33 = pixsize*pixsize;
01304     t41 = 1.0/t7/t6;
01305     t45 = t26*t26;
01306     t47 = t8/t45;
01307     t53 = -t13*t41*fcoll+t15*t8-t5*t15*t41;
01308     t56 = t12/t24;
01309     t57 = t6*t6;
01310     t76 = -xfibre*t15*t8+fcoll*t12*t8;
01311     t93 = gspace*gspace;
01312     t94 = 1.0/t93;
01313 
01314     *y = -t2*t30+0.5*ny;
01315 
01316     /* Check if derivatives expected */
01317     if (dyda == NULL) return;
01318 
01319     /* derivatives for each parameters */
01320     dyda[LMP_NY]    = 0.5;                   /* d(y)/d(ny)      */
01321 
01322     dyda[LMP_PXSIZ] = t2*t28/t33;
01323     dyda[LMP_FCOLL] = -cfact*yfibre*t30+cfact*t5*yfibre*t41*t27*t29+t2*t47*t29*
01324                    (-t53*t15+t56*(2.0*t4/t57*fcoll-2.0*t18*t53)/2.0);
01325     dyda[LMP_CFACT] = -fcoll*yfibre*t30;
01326     dyda[LMP_THETA] = t2*t47*t29*(-t76*t15-t18*t12-t15*t24-t56*t18*t76);
01327     dyda[LMP_ORDER] = t2*t47*t29*(lambda*t10*t15+t56*t18*lambda*t10);
01328     dyda[LMP_SPACE] = t2*t47*t29*(-t9*t94*t15-t56*t18*t9*t94);
01329 
01330 } /* end mrqyoptmod() */
01331 
01404 void
01405 mrqyoptmod2(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
01406             cxdouble dyda[], cxint na)
01407 {
01408 
01409     const cxchar *fctid = "mrqyoptmod2";
01410 
01411     register cxdouble lambda,xfibre,yfibre,pixsize,ny;
01412     /* Optical model parameters */
01413     register cxdouble fcoll,cfact;
01414     /* Grating parameters */
01415     register cxdouble gtheta,gorder,gspace;
01416     /* Slit position parameters */
01417     cxdouble  slitdx,slitdy,slitphi;
01418 
01419     double t1,t102,t103,t11,t112,t117,t118,t12,t123,t13,t136,t14,t141,t145,
01420           t147,t15,t159,t16,t160,t17,t172,t179,t18,t184,t19,t2,t21,t22,t24,
01421           t25,t27,t29,t31,t33,t35,t36,t37,t38,t39,t4,t42,t50,t51,t54,t56,t6,
01422           t62,t65,t66,t68,t7,t85;
01423 
01424     /* check for number of parameters */
01425     if (na != 10) {
01426         cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
01427         return;
01428     }
01429 
01430     *y = 0.0;
01431     if (dyda != NULL) {
01432         dyda[LMP_NY] = dyda[LMP_PXSIZ] = dyda[LMP_FCOLL] = dyda[LMP_CFACT] =
01433                       dyda[LMP_THETA] = dyda[LMP_ORDER] = dyda[LMP_SPACE] =
01434                       dyda[LMP_SOFFX] = dyda[LMP_SOFFY] = dyda[LMP_SPHI] = 0.0;
01435     }
01436 
01437     lambda  = x[LMI_WLEN];
01438     xfibre  = x[LMI_XFIB];
01439     yfibre  = x[LMI_YFIB];
01440 
01441     ny      = a[LMP_NY];
01442     pixsize = a[LMP_PXSIZ];
01443     fcoll   = a[LMP_FCOLL];
01444     cfact   = a[LMP_CFACT];
01445     gtheta  = a[LMP_THETA];
01446     gorder  = a[LMP_ORDER];
01447     gspace  = a[LMP_SPACE];
01448     slitdx  = a[LMP_SOFFX];
01449     slitdy  = a[LMP_SOFFY];
01450     slitphi = a[LMP_SPHI];
01451 
01452     t1 = cfact*fcoll;
01453     t2 = slitphi*slitphi;
01454     t4 = sqrt(1.0-t2);
01455     t6 = yfibre*t4+slitdy;
01456     t7 = t1*t6;
01457     t11 = xfibre*(1.0+slitphi*yfibre)+slitdx;
01458     t12 = t11*t11;
01459     t13 = t6*t6;
01460     t14 = fcoll*fcoll;
01461     t15 = t12+t13+t14;
01462     t16 = sqrt(t15);
01463     t17 = 1/t16;
01464     t18 = lambda*gorder;
01465     t19 = 1/gspace;
01466     t21 = cos(gtheta);
01467     t22 = t11*t21;
01468     t24 = sin(gtheta);
01469     t25 = fcoll*t24;
01470     t27 = -t18*t19+t22*t17+t25*t17;
01471     t29 = 1/t15;
01472     t31 = t27*t27;
01473     t33 = sqrt(1.0-t13*t29-t31);
01474     t35 = -t27*t24+t21*t33;
01475     t36 = 1/t35;
01476     t37 = t17*t36;
01477     t38 = 1/pixsize;
01478     t39 = t37*t38;
01479     t42 = pixsize*pixsize;
01480     t50 = 1/t16/t15;
01481     t51 = t50*t36;
01482     t54 = t35*t35;
01483     t56 = t17/t54;
01484     t62 = -t22*t50*fcoll+t24*t17-t14*t24*t50;
01485     t65 = t21/t33;
01486     t66 = t15*t15;
01487     t68 = t13/t66;
01488     t85 = -t11*t24*t17+fcoll*t21*t17;
01489     t102 = gspace*gspace;
01490     t103 = 1/t102;
01491     t112 = 2.0*t11;
01492     t117 = t21*t17;
01493     t118 = t50*t112;
01494     t123 = t117-t22*t118/2.0-t25*t118/2.0;
01495     t136 = 2.0*t6;
01496     t141 = t50*t136;
01497     t145 = -t22*t141/2.0-t25*t141/2.0;
01498     t147 = t6*t29;
01499     t159 = 1/t4;
01500     t160 = yfibre*t159;
01501     t172 = 2.0*t11*xfibre*yfibre-2.0*t6*yfibre*t159*slitphi;
01502     t179 = t50*t172;
01503     t184 = xfibre*yfibre*t117-t22*t179/2.0-t25*t179/2.0;
01504 
01505     *y = -t7*t39+0.5*ny;
01506 
01507     /* Check if derivatives expected */
01508     if (dyda == NULL) return;
01509 
01510     /* derivatives for each parameters */
01511     dyda[LMP_NY]    = 0.5;                   /* d(y)/d(ny)      */
01512     dyda[LMP_PXSIZ] = t7*t37/t42;            /* d(y)/d(pixsize) */
01513     dyda[LMP_FCOLL] =                        /* d(y)/d(fcoll) */
01514                     -cfact*t6*t39+cfact*t14*t6*t51*t38+
01515                      t7*t56*t38*(-t62*t24+t65*(2.0*t68*fcoll-2.0*t27*t62)/2.0);
01516     dyda[LMP_CFACT] =                        /* d(y)/d(cfact) */
01517                     -fcoll*t6*t39;
01518     dyda[LMP_THETA] =                        /* d(y)/d(gtheta) */
01519                     t7*t56*t38*(-t85*t24-t27*t21-t24*t33-t65*t27*t85);
01520     dyda[LMP_ORDER] =                        /* d(y)/d(gorder) */
01521                     t7*t56*t38*(lambda*t19*t24+t65*t27*lambda*t19);
01522     dyda[LMP_SPACE] =                        /* d(y)/d(gspace) */
01523                     t7*t56*t38*(-t18*t103*t24-t65*t27*t18*t103);
01524     dyda[LMP_SOFFX] =                        /* d(y)/d(slitdx) */
01525                     t7*t51*t38*t112/2.0+t7*t56*t38*(-t123*t24+t65*
01526                      (t68*t112-2.0*t27*t123)/2.0);
01527     dyda[LMP_SOFFY] =                        /* d(y)/d(slitdy) */
01528                     -t1*t39+t7*t51*t38*t136/2.0+t7*t56*t38*(-t145*t24+t65*
01529                      (-2.0*t147+t68*t136-2.0*t27*t145)/2.0);
01530     dyda[LMP_SPHI]  =                        /* d(y)/d(slitphi) */
01531                     t1*t160*slitphi*t17*t36*t38+t7*t51*t38*t172/2.0+
01532                     t7*t56*t38*(-t184*t24+t65*(2.0*t147*t160*slitphi+t68*t172-
01533                                  2.0*t27*t184)/2.0);
01534 
01535 } /* end mrqyoptmod2() */
01536 
01562 void
01563 mrqpsfcos(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
01564           cxdouble dyda[], cxint na)
01565 {
01566 
01567     const cxchar *fctid = "mrqpsfcos";
01568 
01569     cxdouble t1,t10,t13,t14,t15,t16,t2,t26,t3,t4,t5,t6,t7,t8,t9;
01570 
01571     cxdouble amplitude  = a[LMP_AMPL];
01572     cxdouble center     = a[LMP_CENT];
01573     cxdouble background = a[LMP_BACK];
01574     cxdouble width1     = a[LMP_WID1];
01575     cxdouble width2     = a[LMP_WID2];
01576 
01577     /* check for number of parameters */
01578     if (na != 5) {
01579         cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
01580         return;
01581     }
01582 
01583     *y = 0.0;
01584     if (dyda != NULL) {
01585         dyda[LMP_AMPL] = dyda[LMP_CENT] = dyda[LMP_BACK] = dyda[LMP_WID1] =
01586                          dyda[LMP_WID2] = 0.0;
01587     }
01588 
01589     t1 = x[0]-center;
01590     t2 = fabs(t1);
01591     t3 = 1.0/width2;
01592     t4 = t2*t3;
01593     t5 = pow(t4,width1);
01594     t6 = M_PI*t5;
01595     t7 = cos(t6);
01596     t8 = 1.0+t7;
01597     t9 = t8*t8;
01598     t10 = t9*t8;
01599     t13 = amplitude*t9;
01600     t14 = sin(t6);
01601     t15 = t13*t14;
01602     t16 = log(t4);
01603     t26 = (t1>0.0)? 1.0:-1.0;
01604 
01605     if (t2 > width2) {
01606         *y =  background;
01607 
01608         /* Check if derivatives expected */
01609         if (dyda == NULL) return;
01610 
01611         dyda[LMP_AMPL] = dyda[LMP_CENT] = dyda[LMP_BACK] = dyda[LMP_WID1] = 0.0;
01612         dyda[LMP_WID2] = 1.0;
01613     } else {
01614         *y = amplitude*t10/8.0+background;     /* Function value     */
01615 
01616         /* Check if derivatives expected */
01617         if (dyda == NULL)
01618             return;
01619 
01620         dyda[LMP_AMPL] = t10/8.0;                     /* d(y)/d(amplitude)  */
01621                                                /* d(y)/d(center)     */
01622         dyda[LMP_CENT] = 3.0/8.0*t13*t14*M_PI*t5*width1*t26/t2;
01623         dyda[LMP_BACK] = 1.0;                         /* d(y)/d(background) */
01624         dyda[LMP_WID1] = -3.0/8.0*t15*t6*t16;         /* d(y)/d(width1)     */
01625         dyda[LMP_WID2] = 3.0/8.0*t15*t6*width1*t3;    /* d(y)/d(width2)     */
01626     }
01627 
01628     return;
01629 
01630 } /* end mrqpsfcos() */
01631 
01657 void
01658 mrqpsfexp(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
01659           cxdouble dyda[], cxint na)
01660 {
01661 
01662     const cxchar *fctid = "mrqpsfexp";
01663 
01664     cxdouble t1,t2,t3,t4,t6,t8,t10,t15,t18;
01665 
01666     cxdouble amplitude  = a[LMP_AMPL];
01667     cxdouble center     = a[LMP_CENT];
01668     cxdouble background = a[LMP_BACK];
01669     cxdouble width1     = a[LMP_WID1];
01670     cxdouble width2     = a[LMP_WID2];
01671 
01672     /* check for number of parameters */
01673     if (na != 5) {
01674         cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
01675         return;
01676     }
01677 
01678     *y = 0.0;
01679     if (dyda != NULL) {
01680         dyda[LMP_AMPL] = dyda[LMP_CENT] = dyda[LMP_BACK] = dyda[LMP_WID1] =
01681                          dyda[LMP_WID2] = 0.0;
01682     }
01683 
01684     t1 = x[0]-center;
01685 
01686     if (t1 > 0.0) {
01687         t2 = t1;
01688         t10 = 1.0;
01689     } else {
01690         t2 = -t1;
01691         t10 = -1.0;
01692     }
01693 
01694     t3 = pow(t2,width2);
01695     t4 = 1.0/width1;
01696     t6 = exp(-t3*t4);
01697     t8 = amplitude*t3;
01698     t15 = width1*width1;
01699     t18 = log(t2);
01700 
01701     *y = amplitude*t6+background;
01702 
01703     /* Check if derivatives expected */
01704     if (dyda == NULL)
01705         return;
01706 
01707     dyda[LMP_AMPL] = t6;                      /* d(y)/d(amplitude)  */
01708     dyda[LMP_CENT] = t8*width2*t10/t2*t4*t6;  /* d(y)/d(center)     */
01709 
01710     if (isnan(dyda[LMP_CENT]))
01711         dyda[LMP_CENT] = 0.;
01712 
01713     dyda[LMP_BACK] = 1.0;                     /* d(y)/d(background) */
01714     dyda[LMP_WID1] = t8/t15*t6;               /* d(y)/d(width1)     */
01715     dyda[LMP_WID2] = -t8*t18*t4*t6;           /* d(y)/d(width2)     */
01716 
01717     if (isnan(dyda[LMP_WID2]))
01718         dyda[LMP_WID2] = 0.;
01719 
01720     if (r != NULL) {
01721         register cxint k;
01722 
01723         k = LMP_AMPL << 1;
01724         if (r[k+1] > 0) {
01725             dyda[LMP_AMPL] *= mrqdydaweight(a[LMP_AMPL],r[k],r[k+1]);
01726         }
01727         k = LMP_CENT << 1;
01728         if (r[k+1] > 0) {
01729             dyda[LMP_CENT] *= mrqdydaweight(a[LMP_CENT],r[k],r[k+1]);
01730         }
01731         k = LMP_WID1 << 1;
01732         if (r[k+1] > 0) {
01733             dyda[LMP_WID1] *= mrqdydaweight(a[LMP_WID1],r[k],r[k+1]);
01734         }
01735         k = LMP_WID2 << 1;
01736         if (r[k+1] > 0) {
01737             dyda[LMP_WID2] *= mrqdydaweight(a[LMP_WID2],r[k],r[k+1]);
01738         }
01739     }
01740 
01741     return;
01742 
01743 } /* end mrqpsfexp() */
01744 
01770 void
01771 mrqpsfexp2(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
01772            cxdouble dyda[], cxint na)
01773 {
01774 
01775     const cxchar *fctid = "mrqpsfexp2";
01776 
01777     cxdouble t1,t2,t3,t4,t5,t6,t8,t10,t16;
01778 
01779     cxdouble amplitude  = a[LMP_AMPL];
01780     cxdouble center     = a[LMP_CENT];
01781     cxdouble background = a[LMP_BACK];
01782     cxdouble width1     = a[LMP_WID1];
01783     cxdouble width2     = a[LMP_WID2];
01784 
01785     /* check for number of parameters */
01786     if (na != 5) {
01787         cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
01788         return;
01789     }
01790 
01791     *y = 0.0;
01792     if (dyda != NULL) {
01793         dyda[LMP_AMPL] = dyda[LMP_CENT] = dyda[LMP_BACK] = dyda[LMP_WID1] =
01794                          dyda[LMP_WID2] = 0.0;
01795     }
01796 
01797     t1 = x[0]-center;
01798 
01799     if (t1 > 0.0) {
01800         t2 = t1;
01801         t10 = 1.0;
01802     } else {
01803         t2 = -t1;
01804         t10 = -1.0;
01805     }
01806 
01807     t3 = 1.0/width1;
01808     t4 = t2*t3;
01809     t5 = pow(t4,width2);
01810     t6 = exp(-t5);
01811     t8 = amplitude*t5;
01812     t16 = log(t4);
01813 
01814     *y = amplitude*t6+background;
01815 
01816     /* Check if derivatives expected */
01817     if (dyda == NULL)
01818         return;
01819 
01820     dyda[LMP_AMPL] = t6;                         /* d(y)/d(amplitude)  */
01821     dyda[LMP_CENT] = t8*width2*t10/t2*t6;        /* d(y)/d(center)     */
01822 
01823     if (isnan(dyda[LMP_CENT]))
01824         dyda[LMP_CENT] = 0.0;
01825 
01826     dyda[LMP_BACK] = 1.0;                        /* d(y)/d(background) */
01827     dyda[LMP_WID1] = t8*width2*t3*t6;            /* d(y)/d(width1)     */
01828     dyda[LMP_WID2] = -t8*t16*t6;                 /* d(y)/d(width2)     */
01829 
01830     if (isnan(dyda[LMP_WID2]))
01831         dyda[LMP_WID2] = 0.0;
01832 
01833     if (r != NULL) {
01834         register cxint k;
01835 
01836         k = LMP_AMPL << 1;
01837         if (r[k+1] > 0) {
01838             dyda[LMP_AMPL] *= mrqdydaweight(a[LMP_AMPL],r[k],r[k+1]);
01839         }
01840         k = LMP_CENT << 1;
01841         if (r[k+1] > 0) {
01842             dyda[LMP_CENT] *= mrqdydaweight(a[LMP_CENT],r[k],r[k+1]);
01843         }
01844         k = LMP_WID1 << 1;
01845         if (r[k+1] > 0) {
01846             dyda[LMP_WID1] *= mrqdydaweight(a[LMP_WID1],r[k],r[k+1]);
01847         }
01848         k = LMP_WID2 << 1;
01849         if (r[k+1] > 0) {
01850             dyda[LMP_WID2] *= mrqdydaweight(a[LMP_WID2],r[k],r[k+1]);
01851         }
01852     }
01853 
01854     return;
01855 
01856 } /* end mrqpsfexp2() */
01857 
01878 void
01879 mrqlocywarp(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
01880             cxdouble dyda[], cxint na)
01881 {
01882 
01883     const cxchar *fctid = "mrqlocywarp";
01884 
01885     cxdouble          xccd, nx, startx;
01886     register cxint    i,ncoef;
01887     cxdouble          Tx, Ty, cx, Ky, tt, xx;
01888     cpl_matrix       *mBase = NULL, *mX = NULL;
01889     register cxdouble fxx = 0.0, f1xx = 0.0, f2xx = 0.0, z1;
01890     cxdouble         *pd_x = NULL, *pd_mbase = NULL;
01891 
01892     /* check for number of parameters */
01893     if (na != 5) {
01894         cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
01895         return;
01896     }
01897 
01898     *y = 0.0;
01899     if (dyda != NULL) {
01900         dyda[LMP_TX] = dyda[LMP_TY] = dyda[LMP_CX] = dyda[LMP_KY] =
01901                        dyda[LMP_TT] = 0.0;
01902     }
01903 
01904     xccd   = x[LMI_XCCD];        /* pixel abcissa */
01905     nx     = x[LMI_NX];          /* number of pixels in dispersion dir. */
01906     startx = x[LMI_STRX];        /* 1st pixel position */
01907     ncoef  = (cxint) x[LMI_NCOF];  /* number of chebyshev coef */
01908 
01909     Tx = a[LMP_TX];   /* Translation X */
01910     Ty = a[LMP_TY];   /* Translation Y */
01911     cx = a[LMP_CX];   /* Scaling X: cx = 1/Kx */
01912     Ky = a[LMP_KY];   /* Scaling Y */
01913     tt = a[LMP_TT];   /* Rotation theta: tt = tan(theta) */
01914 
01915     xx = cx*(xccd-Tx);
01916 
01917     mX = cpl_matrix_new(1,1);
01918     pd_x = cpl_matrix_get_data(mX);
01919     pd_x[0] = xx;
01920 
01921     mBase = giraffe_chebyshev_base1d(startx, nx, ncoef, mX);
01922 
01923     pd_mbase = cpl_matrix_get_data(mBase);
01924 
01925     for (i = 0; i < ncoef; i++)
01926         fxx += pd_mbase[i] * x[i+4];
01927 
01928     if (ncoef > 1) {
01929         for (i = 0; i < (ncoef - 1); i++)
01930             f1xx += pd_mbase[i] * (i+1)*x[i+5];
01931     }
01932 
01933     if (ncoef > 2) {
01934         for (i = 0; i < (ncoef - 2); i++)
01935             f2xx += pd_mbase[i] * (i+2)*x[i+6];
01936     }
01937 
01938     if (mX!=NULL) { cpl_matrix_delete(mX); mX = NULL; pd_x = NULL; }
01939     if (mBase!=NULL) { cpl_matrix_delete(mBase); mBase = NULL; pd_mbase = NULL; }
01940 
01941     z1 = 1.0 - tt*tt + tt*f1xx;
01942     *y = Ky*(fxx-tt*xx)/z1 + Ty;
01943 
01944     /* Check if derivatives expected */
01945     if (dyda == NULL)
01946         return;
01947 
01948     dyda[LMP_TX] =                          /* d(y)/d(Tx)  */
01949               (cx*Ky/z1)*((tt-f1xx) + tt*f2xx*(fxx-tt*xx)/z1);
01950 
01951     dyda[LMP_TY] = 1.0;                     /* d(y)/d(Ty)     */
01952 
01953     dyda[LMP_CX] =                          /* d(y)/d(cx)     */
01954               (Ky*(xccd-Tx)/z1)*((f1xx-tt) - tt*f2xx*(fxx-tt*xx)/z1);
01955 
01956     dyda[LMP_KY] = (fxx-tt*xx)/z1;              /* d(y)/d(Ky)     */
01957     dyda[LMP_TT] =                          /* d(y)/d(tt) */
01958               (Ky/(z1*z1))*(-xx*(1.+tt*tt)+2.*tt*fxx-fxx*f1xx);
01959 
01960     if (r != NULL) {
01961         register cxint k;
01962 
01963         k = LMP_TX << 1;
01964         if (r[k+1] > 0) {
01965             dyda[LMP_TX] *= mrqdydaweight(a[LMP_TX],r[k],r[k+1]);
01966         }
01967         k = LMP_CX << 1;
01968         if (r[k+1] > 0) {
01969             dyda[LMP_CX] *= mrqdydaweight(a[LMP_CX],r[k],r[k+1]);
01970         }
01971         k = LMP_KY << 1;
01972         if (r[k+1] > 0) {
01973             dyda[LMP_KY] *= mrqdydaweight(a[LMP_KY],r[k],r[k+1]);
01974         }
01975         k = LMP_TT << 1;
01976         if (r[k+1] > 0) {
01977             dyda[LMP_TT] *= mrqdydaweight(a[LMP_TT],r[k],r[k+1]);
01978         }
01979     }
01980 
01981     return;
01982 
01983 } /* end mrqlocywarp() */
01984 
02005 void
02006 mrqxoptmodGS(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
02007              cxdouble dyda[], cxint na)
02008 {
02009 
02010     const cxchar *fctid = "mrqxoptmodGS";
02011 
02012     register cxdouble lambda,xfibre,yfibre,pixsize,nx;
02013     /* Optical model parameters */
02014     register cxdouble fcoll,cfact;
02015     /* Grating parameters */
02016     register cxdouble gtheta,gorder,gspace;
02017 
02018     register cxdouble t1,t10,t109,t11,t110,t114,t12,t14,t17,t18,t2,t21,t23,t24;
02019     register cxdouble t25,t26,t27,t28,t29,t3,t30,t31,t35,t40,t43,t49,t5,t51,t52;
02020     register cxdouble t53,t59,t6,t66,t67,t69,t7,t71,t8,t82,t84,t9,t95,t98;
02021 
02022     /* check for number of parameters */
02023     if (na != 7) {
02024         cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
02025         return;
02026     }
02027 
02028     *y = 0.0;
02029     if (dyda != NULL) {
02030         dyda[LMP_NX] = dyda[LMP_PXSIZ] = dyda[LMP_FCOLL] = dyda[LMP_CFACT] =
02031                       dyda[LMP_THETA] = dyda[LMP_ORDER] = dyda[LMP_SPACE] = 0.0;
02032     }
02033 
02034     lambda  = x[LMI_WLEN];    /* wavelength [mm]              */
02035     xfibre  = x[LMI_XFIB];    /* Y fibre [mm]                 */
02036     yfibre  = x[LMI_YFIB];    /* Y fibre [mm]                 */
02037 
02038     nx      = a[LMP_NX];      /* CCD  size in X [pixels]      */
02039     pixsize = a[LMP_PXSIZ];   /* CCD pixel size [mm]          */
02040     fcoll   = a[LMP_FCOLL];   /* collimator focal length [mm] */
02041     cfact   = a[LMP_CFACT];   /* camera magnification factor  */
02042     gtheta  = a[LMP_THETA];   /* grating angle [radian]       */
02043     gorder  = a[LMP_ORDER];   /* grating diffraction order    */
02044     gspace  = a[LMP_SPACE];   /* grating groove spacing [mm]  */
02045 
02046     t1 = cfact*fcoll;
02047     t2 = lambda*gorder;
02048     t3 = 1.0/gspace;
02049     t5 = cos(gtheta);
02050     t6 = xfibre*t5;
02051     t7 = xfibre*xfibre;
02052     t8 = yfibre*yfibre;
02053     t9 = fcoll*fcoll;
02054     t10 = t7+t8+t9;
02055     t11 = sqrt(t10);
02056     t12 = 1.0/t11;
02057     t14 = sin(gtheta);
02058     t17 = -t2*t3+t12*t6+fcoll*t14*t12;
02059     t18 = t17*t5;
02060     t21 = t17*t17;
02061     t23 = sqrt(1.0-t8/t10-t21);
02062     t24 = t14*t23;
02063     t25 = t18+t24;
02064     t26 = t17*t14;
02065     t27 = t5*t23;
02066     t28 = -t26+t27;
02067     t29 = 1.0/t28;
02068     t30 = t25*t29;
02069     t31 = 1.0/pixsize;
02070     t35 = pixsize*pixsize;
02071     t40 = t29*t31;
02072     t43 = 1.0/t11/t10;
02073     t49 = -t6*t43*fcoll+t14*t12-t9*t14*t43;
02074     t51 = 1.0/t23;
02075     t52 = t14*t51;
02076     t53 = t10*t10;
02077     t59 = 2.0*t8/t53*fcoll-2.0*t17*t49;
02078     t66 = t1*t25;
02079     t67 = t28*t28;
02080     t69 = 1.0/t67*t31;
02081     t71 = t5*t51;
02082     t82 = -xfibre*t14*t12+fcoll*t5*t12;
02083     t84 = t17*t82;
02084     t95 = lambda*t3;
02085     t98 = t17*lambda*t3;
02086     t109 = gspace*gspace;
02087     t110 = 1.0/t109;
02088     t114 = t2*t110;
02089 
02090     /* takes care of model direction */
02091     if (nx < 0.0)
02092         *y = t1*t30*t31-0.5*nx;
02093     else
02094         *y = -t1*t30*t31+0.5*nx;
02095 
02096     /* Check if derivatives expected */
02097     if (dyda == NULL)
02098         return;
02099 
02100     /* derivatives for each parameters */
02101     dyda[LMP_NX]    = 0.5;                   /* d(y)/d(nx)      */
02102     dyda[LMP_PXSIZ] = -t1*t30/t35;           /* d(y)/d(pixsize) */
02103     dyda[LMP_FCOLL] =                        /* d(y)/d(fcoll) */
02104                     cfact*t25*t40+t1*(t49*t5+t52*t59/2.0)*t29*t31-
02105                     t66*t69*(-t49*t14+t71*t59/2.0);
02106     dyda[LMP_CFACT] =                        /* d(y)/d(cfact) */
02107                     fcoll*t25*t40;
02108     dyda[LMP_THETA] =                        /* d(y)/d(gtheta) */
02109                     t1*(t82*t5-t26+t27-t52*t84)*t29*t31-
02110                     t66*t69*(-t82*t14-t18-t24-t71*t84);
02111     dyda[LMP_ORDER] =                        /* d(y)/d(gorder) */
02112                     t1*(-t95*t5+t52*t98)*t29*t31-t66*t69*(t95*t14+t71*t98);
02113     dyda[LMP_SPACE] =                        /* d(y)/d(gspace) */
02114                     t1*(t2*t110*t5-t52*t17*t114)*t29*t31-
02115                     t66*t69*(-t2*t110*t14-t71*t17*t114);
02116 
02117     if (nx > 0.0) {
02118         dyda[LMP_NX]    = -dyda[LMP_NX];
02119         dyda[LMP_PXSIZ] = -dyda[LMP_PXSIZ];
02120         dyda[LMP_FCOLL] = -dyda[LMP_FCOLL];
02121         dyda[LMP_CFACT] = -dyda[LMP_CFACT];
02122         dyda[LMP_THETA] = -dyda[LMP_THETA];
02123         dyda[LMP_ORDER] = -dyda[LMP_ORDER];
02124         dyda[LMP_SPACE] = -dyda[LMP_SPACE];
02125     }
02126 
02127     if (r != NULL) {
02128         register cxint k;
02129 
02130         k = LMP_PXSIZ << 1;
02131         if (r[k+1] > 0) {
02132             dyda[LMP_PXSIZ] *= mrqdydaweight(a[LMP_PXSIZ],r[k],r[k+1]);
02133         }
02134         k = LMP_FCOLL << 1;
02135         if (r[k+1] > 0) {
02136             dyda[LMP_FCOLL] *= mrqdydaweight(a[LMP_FCOLL],r[k],r[k+1]);
02137         }
02138         k = LMP_CFACT << 1;
02139         if (r[k+1] > 0) {
02140             dyda[LMP_CFACT] *= mrqdydaweight(a[LMP_CFACT],r[k],r[k+1]);
02141         }
02142         k = LMP_THETA << 1;
02143         if (r[k+1] > 0) {
02144             dyda[LMP_THETA] *=  mrqdydaweight(a[LMP_THETA],r[k],r[k+1]);
02145         }
02146         k = LMP_ORDER << 1;
02147         if (r[k+1] > 0) {
02148             dyda[LMP_ORDER] *=  mrqdydaweight(a[LMP_ORDER],r[k],r[k+1]);
02149         }
02150         k = LMP_SPACE << 1;
02151         if (r[k+1] > 0) {
02152             dyda[LMP_SPACE] *=  mrqdydaweight(a[LMP_SPACE],r[k],r[k+1]);
02153         }
02154     }
02155 
02156 } /* end mrqxoptmodGS() */
02157 
02181 void
02182 mrqtest(cxdouble x[], cxdouble a[], cxdouble r[], cxdouble *y,
02183            cxdouble dyda[], cxint na)
02184 {
02185 
02186     const cxchar *fctid = "mrqtest";
02187 
02188     cxdouble a1 = a[0];
02189     cxdouble b1 = a[1];
02190 
02191     /* check for number of parameters */
02192     if (na != 2) {
02193         cpl_error_set(fctid, CPL_ERROR_ILLEGAL_INPUT);
02194         return;
02195     }
02196 
02197     *y = 0.0;
02198     *y = a1 * x[0] + b1;
02199 
02200     /* Check if derivatives expected */
02201     if (dyda == NULL)
02202         return;
02203 
02204     dyda[0] = 0.0;
02205     dyda[1] = 0.0;
02206 
02207     return;
02208 
02209 } /* end mrqtest() */
02210 

This file is part of the GIRAFFE Pipeline Reference Manual 2.8.8.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Fri Mar 4 10:50:27 2011 by doxygen 1.6.3 written by Dimitri van Heesch, © 1997-2004