scalar XI_CylShell(scalar x, sasfit_param * param)
{
	scalar xR, xH, XIres;

	SASFIT_ASSERT_PTR(param);

	xR = x*R*sin(ALPHA);
	xH = x*H*cos(ALPHA)/2.;

	if ((R+H) == 0) return 0;

	if (xR == 0)
	{
		XIres = R/(R+H) * cos(xH);
	} else
	{
		XIres = R/(R+H) * 2.*gsl_sf_bessel_J1(xR)/xR * cos(xH);
	}

	if (xH == 0)
	{
		XIres = XIres + H/(R+H) * gsl_sf_bessel_J0(xR);
	} else
	{
		XIres = XIres + H/(R+H) * gsl_sf_bessel_J0(xR)*sin(xH)/xH;
	}
	return (2.* M_PI*R*R + 2.0*M_PI*R*H)*XIres;
}
scalar sasfit_ff_very_long_cyl_w_gauss_chains(scalar q, sasfit_param * param)
{
	scalar Pp, Pcs, Fc, w, Scc, Ssc, J1x_x;
	scalar R, Rg, d, Nagg, H, r_core, r_chains;

	SASFIT_ASSERT_PTR(param);

	sasfit_get_param(param, 7, &R, &Rg, &d, &Nagg, &H, &r_core, &r_chains);

	SASFIT_CHECK_COND1((H <= 0.0), param, "H(%lg) < 0", H);

	if (q == 0.0) return Pp = 1.0;

	Pp =(2.0*gsl_sf_Si(q*H)/(q*H)-pow(sin(q*H*0.5)/(q*H*0.5),2.0));

	w = sasfit_rwbrush_w(q, Rg);
	Fc = sasfit_gauss_fc(q, Rg);

	Scc = w*w*pow(gsl_sf_bessel_J0(q*(R+d*Rg)),2.);
	if (q*R == 0)
	{
		J1x_x = 0.5;
	} else
	{
		J1x_x = gsl_sf_bessel_J1(q*R)/(q*R);
	}
	Ssc = w*2.*J1x_x*gsl_sf_bessel_J0(q*(R+d*Rg));

	Pcs = pow(r_core*2.*J1x_x,2.)
		+ Nagg*(Nagg-1.)*pow(r_chains,2.0)*Scc*((Nagg < 1) ?  0 : 1)
		+ 2.*Nagg*r_chains*r_core*Ssc;
//  sasfit_out("Cyl_Fc:%lf \n",pow(r_chains,2.0)*Nagg*Fc);
	return Pp*Pcs
		+ pow(r_chains,2.0)*Nagg*Fc;
}
Example #3
0
scalar XI(scalar Q, scalar R, scalar H, scalar alpha)
{
        scalar xR, xH, XIres;

        xR = Q*R*sin(alpha);
        xH = Q*H*cos(alpha)/2.0;

        if ((R+H) == 0.0) return 0.0;

        if (xR == 0.0)
	{
                XIres = R/(R+H) * cos(xH);
        } else {
                XIres = R/(R+H) * 2.0*gsl_sf_bessel_J1(xR)/xR * cos(xH);
        }

        if (xH == 0.0)
	{
                XIres = XIres + H/(R+H) * gsl_sf_bessel_J0(xR);
        } else {
                XIres = XIres + H/(R+H) * gsl_sf_bessel_J0(xR)*sin(xH)/xH;
        }

        return XIres;
}
const Real GreensFunction2DAbs::p_int_theta_first(const Real r,
                                                  const Real theta,
                                                  const Real t) const
{
    const Real r_0(this->getr0());
    const Real a(this->geta());
    const Real minusDt(-1e0 * this->getD() * t);

    const Integer num_term_use(100);
    const Real threshold(CUTOFF);

    Real sum(0e0);
    Real term(0e0);
    Real term1(0e0);
    Real term2(0e0);
    Real term3(0e0);

    Real a_alpha_n(0e0);
    Real alpha_n(0e0);
    Real J0_r_alpha_n(0e0);
    Real J0_r0_alpha_n(0e0);
    Real J1_a_alpha_n(0e0);

    Integer n(1);
    for(; n < num_term_use; ++n)
    {
        a_alpha_n = gsl_sf_bessel_zero_J0(n);
        alpha_n = a_alpha_n / a;
        J0_r_alpha_n  = gsl_sf_bessel_J0(r * alpha_n);
        J0_r0_alpha_n = gsl_sf_bessel_J0(r_0 * alpha_n);
        J1_a_alpha_n  = gsl_sf_bessel_J1(a_alpha_n);

        term1 = std::exp(alpha_n * alpha_n * minusDt);
        term2 = J0_r_alpha_n * J0_r0_alpha_n;
        term3 = J1_a_alpha_n * J1_a_alpha_n;

        term = term1 * term2 / term3;
        sum += term;

//             std::cout << "sum " << sum << ", term" << term << std::endl;

        if(fabs(term/sum) < threshold)
        {
//                 std::cout << "normal exit. n = " << n << " first term" << std::endl;
            break;
        }
    }
    if(n == num_term_use)
        std::cout << "warning: use term over num_term_use" << std::endl;

//         return (sum / (M_PI * a * a));
    return (theta * sum / (M_PI * a * a));
}
Example #5
0
static double lenscfp_integrand(double ell, void *p)
{
  lensCorrFunc lcf = (lensCorrFunc) p;
  
  return ell/2.0/M_PI*lens_power_spectrum(ell,lcf->lps)*gsl_sf_bessel_J0(ell*lcf->theta/60.0/180.0*M_PI);
  
}
Example #6
0
/*
 * Class:     GSLOps_GSLOps
 * Method:    gslsfbesselJ0
 * Signature: (D)D
 */
JNIEXPORT jdouble JNICALL Java_GSLJava_GSLOps_gslSfBesselJ0
  (JNIEnv *env, jobject obj, jdouble  x)
 {

  double y = gsl_sf_bessel_J0(x);	
 return y;
}
Example #7
0
/*
 * form factor of a spherical Cylinder with radius R, height L and scattering
 * length density eta
 */
scalar sasfit_ff_long_cyl(scalar q, sasfit_param * param)
{
	scalar mu, sigma, pi_mu, G1, G2, I_sp, Sum, V, omega;
	scalar R, L, eta;

	SASFIT_ASSERT_PTR(param);

	sasfit_get_param(param, 4, &R, &L, EMPTY, &eta);

	SASFIT_CHECK_COND1((q < 0.0), param, "q(%lg) < 0",q);
	SASFIT_CHECK_COND1((R < 0.0), param, "R(%lg) < 0",R);
	SASFIT_CHECK_COND1((L < 0.0), param, "L(%lg) < 0",L);

	if (R == 0.0) return 0.0;
	if (L == 0.0) return 0.0;

	mu = L*q;
	sigma = 2.0*R*q;
	V = M_PI*R*R*L;

	if (R==0 || L==0) return 0;
	if (q==0) return V*V*eta*eta;

	pi_mu = gsl_sf_Si(mu)+cos(mu)/mu+sin(mu)/mu/mu;
	G1 = 2.0/(0.5*sigma) *  gsl_sf_bessel_J1(0.5*sigma);
	G2 = 8.0/sigma/sigma * gsl_sf_bessel_Jn(2,sigma);
//	I_sp = 3.0 * (sin(sigma*0.5)-0.5*sigma*cos(0.5*sigma)) / pow(sigma/2.0,3);
//	I_sp = I_sp*I_sp;
	omega = 8/gsl_pow_2(sigma)*(3*gsl_sf_bessel_Jn(2,sigma)+gsl_sf_bessel_J0(sigma)-1);

//	Sum = 2.0/mu * (pi_mu*G1*G1 - 1.0/mu*(2.0*G2-I_sp) - sin(mu)/mu/mu);
	Sum = 2.0/mu * (pi_mu*G1*G1 - omega/mu - sin(mu)/mu/mu);

	return eta*eta *V*V* Sum;
}
Example #8
0
void FELNumerical::initParams()
{
    au     = 0.934*bfield*100*lambdau/sqrt(2);      // normalized undulator parameter, calculated from the input parameters read from namelist file
    ku     = 2*PI/lambdau;                          // undulator wave number
    omegas = 2*PI*C0/lambdas;                       // radiation angular frequency
    gammar = sqrt(lambdau*(1+au*au)/2.0/lambdas);   // resonant gamma
    sigmax = sqrt(avgbeta*emitn/gamma0);
    j0     = current/(PI*sigmax*sigmax);            // transverse current density
    coef1  = E0/M0/C0/C0;
    coef2  = mu0*C0*C0/omegas;
    coef3  = mu0*C0/4.0/gammar;
    ndelz  = lambdau*deltz;                         // integration step size, [m]
    totalIntSteps = (int) (num/deltz);              // total integration steps
    K0Arr   = new double [totalIntSteps];
    JJArr   = new double [totalIntSteps];
    zposArr = new double [totalIntSteps];
    bfArr   = new double [totalIntSteps];
    ExArr   = new efield [totalIntSteps];
    EyArr   = new efield [totalIntSteps];

    double btmp;
    for(unsigned int i = 0; i < totalIntSteps; i++)
    {
        K0Arr[i] = au*sqrt(2);
        btmp = K0Arr[i]*K0Arr[i]/(4.0+2.0*K0Arr[i]*K0Arr[i]);
        JJArr[i] = gsl_sf_bessel_J0(btmp) - gsl_sf_bessel_J1(btmp);
    }
}
scalar sasfit_Rod_Exp_ave_core(scalar x, sasfit_param * param)
{
    scalar Q, u, Rc, alpha, etashin, etashout, etaexp, t, arg;

    SASFIT_ASSERT_PTR(param);

    sasfit_get_param(param, 4, &Rc, &etashout, &t, &alpha);
    Q		= param->p[MAXPAR-1];

    u=Q*x;
    etashin = 1.0;
    arg = (x-Rc)/t;

    if (x<Rc)	return 0.0;
    if (x>Rc+t)	return 0.0;

    if (alpha < 0)
    {
        etaexp = etashin+(etashout-etashin)*arg*exp((1-arg)*alpha);
    } else
    {
        etaexp = etashout+(etashin-etashout)*(1-arg)*exp(-arg*alpha);
    }

    if (u==0.0) return etaexp * 2.0 * M_PI * x;

    return etaexp * 2.0 * M_PI * x* gsl_sf_bessel_J0(u);
}
Example #10
0
int
main (void)
{
  double x = 5.0;
  double y = gsl_sf_bessel_J0 (x);
  printf ("J0(%g) = %.18e\n", x, y);
  return 0;
}
scalar sasfit_ff_worm_w_gauss_chains(scalar q, sasfit_param * param)
{
	scalar Pp, Pcs, Fc,w, Scc, Ssc, J1x_x;
	scalar R, Rg, d, Nagg, l, L, r_core, r_chains;
	sasfit_param subParam;

	SASFIT_ASSERT_PTR(param);

	sasfit_get_param(param, 8, &R, &Rg, &d, &Nagg, &l, &L, &r_core, &r_chains);

	SASFIT_CHECK_COND1((L <= 0.0), param, "L(%lg) < 0", L);


	//Pp = KholodenkoWorm(interp,q,0.0,l,L,error);
	if (q==0) {
		Pp = 1;
		w = 1;
		Fc = 1;
	} else {
		sasfit_init_param(&subParam);
		subParam.p[0] = 0.0;
		subParam.p[1] = l;
		subParam.p[2] = L;
		Pp = sasfit_ff_KholodenkoWorm(q, &subParam);
		SASFIT_CHECK_SUB_ERR(param, subParam);

		w = sasfit_rwbrush_w(q, Rg);
		Fc = sasfit_gauss_fc(q, Rg);
	}
	Scc = w*w*pow(gsl_sf_bessel_J0(q*(R+d*Rg)),2.);
	if (q*R == 0)
	{
		J1x_x =0.5;
	} else
	{
		J1x_x = gsl_sf_bessel_J1(q*R)/(q*R);
	}
	Ssc = w*2.*J1x_x*gsl_sf_bessel_J0(q*(R+d*Rg));

	Pcs = pow(r_core*2.*J1x_x,2.)
		+ Nagg*(Nagg-1.)*pow(r_chains,2.0)*Scc*((Nagg < 1) ?  0 : 1)
		+ 2.*Nagg*r_chains*r_core*Ssc;

	return Pp*Pcs
		+ pow(r_chains,2.0)*Nagg*Fc;
}
int main(int argc, char** argv) {

    double x = 5.0;
    double y = gsl_sf_bessel_J0(x);
    printf("J0(%5.5le) = %.18le\n", x, y);
    
    return (EXIT_SUCCESS);
}
const Real GreensFunction2DAbs::p_int_r(const Real r, const Real t) const
{
    //speed of convergence is too slow
    if(r == 0e0) return 0e0;

    const Real r_0(this->getr0());
    const Real a(this->geta());
    const Real Dt(this->getD() * t);
    const Integer num_term_use(100);
    const Real threshold(CUTOFF);

    Real sum(0e0);
    Real term(0e0);
    Real term1(0e0);
    Real term2(0e0);
    Real term3(0e0);

    Real a_alpha_n(0e0);
    Real alpha_n(0e0);
    Real J0_r0_alpha_n(0e0);
    Real J1_r_alpha_n(0e0);
    Real J1_a_alpha_n(0e0);

    Integer n(1);
    for(; n < num_term_use; ++n)
    {
        a_alpha_n = gsl_sf_bessel_zero_J0(n);
        alpha_n = a_alpha_n / a;
        J0_r0_alpha_n = gsl_sf_bessel_J0(r_0 * alpha_n);
        J1_r_alpha_n  = gsl_sf_bessel_J1(r * alpha_n);
        J1_a_alpha_n  = gsl_sf_bessel_J1(a_alpha_n);

        term1 = std::exp(-1e0 * alpha_n * alpha_n * Dt);
        term2 = r * J1_r_alpha_n * J0_r0_alpha_n;
        term3 = (alpha_n * J1_a_alpha_n * J1_a_alpha_n);

        term = term1 * term2 / term3;
        sum += term;

//             std::cout << "sum " << sum << ", term" << term << std::endl;

        if(fabs(term/sum) < threshold)
        {
//                 std::cout << "normal exit. " << n << std::endl;
            break;
        }
    }
    if(n == num_term_use)
        std::cout << "warning: use term over num_term_use" << std::endl;

    return (2e0 * sum / (a*a));
}
Example #14
0
static double f_kx_zero(double w, void *p)
{
  struct parameters *params = (struct parameters *)p;
  double  kr = (params->kr);
  double  x0 = (params->x0);
  double result;
  
  if (w == 0.0) w = zero;
  
  result = 2.0*M_PI*w*gsl_sf_bessel_J0(kr*w)*
    (log(x0 + sqrt(w*w + x0*x0)) - log(-x0 + sqrt(w*w + x0*x0)));
  return result;
}
Example #15
0
double BesselFunction::theta_callback(double theta, void *params) 
{
    IntegrationInfo *i = static_cast<IntegrationInfo*>(params);

    double cos_theta = cos(theta);
    std::complex<double> rv = std::complex<double>( 
        sqrt(cos_theta) * sin(theta) * i->wavenumber.value() * gsl_sf_bessel_J0( sin(theta) * i->bessel_factor ), 0.0 )
        * exp( std::complex<double>(0, cos_theta * i->z_offset ) );

    if ( imag_part )
        return imag(rv);
    else
        return real(rv);
}
Example #16
0
/*
  input: alpha:   k*dx (wavenumber * grid spacing)
         points:  points to evaluate bessles at
         npoints: number of points in points
	 M:       highest order bessel function to use

  output: return value: matrix A where
          A[i][j] = J_{j'}(alpha*r_i) /  * f(j'*\theta_i) where
	    0 <= i < npoints
	    0 <= j <= M
	    j' = j     if j <= M
	         j - M if j > M
	    f = 1   if j = 0
	        sin if 1 <= j <= M
	        cos if j > M
 */
gsl_matrix *bessel_matrix(double alpha, point *points, int npoints, int M, double r_typical) {
  gsl_matrix *out = gsl_matrix_alloc(npoints, 2*M+1);
  int i,j;
  double x, y, r, theta;

  double column_scale;

  for (i = 0 ; i < npoints ; i++) { // loop over rows
    x = points[i].x;
    y = points[i].y;
    r = R(x, y);
    theta = THETA(x,y);

    // loops over columns
    gsl_matrix_set(out, i, 0, gsl_sf_bessel_J0(alpha * r) / gsl_sf_bessel_J0(alpha * r_typical));
    for (j = 1 ; j <= M ; j++) {
      column_scale = gsl_sf_bessel_Jn(j, alpha * r_typical);
      gsl_matrix_set(out, i, j, gsl_sf_bessel_Jn(j, alpha * r) * sin(j * theta) / column_scale);
      gsl_matrix_set(out, i, j+M, gsl_sf_bessel_Jn(j, alpha * r) * cos(j * theta) / column_scale);
    }
  }

  return out;
}
const Real GreensFunction2DAbs::p_survival(const Real t) const
{
    // when t == 0.0, return value become eventually 1.0,
    // but the speed of convergence is too slow.
    if(t == 0.0) return 1.0;

    const Real r_0(this->getr0());
    const Real a(this->geta());
    const Real Dt(this->getD() * t);
    const Integer num_term_use(100);
    const Real threshold(CUTOFF);

    Real sum(0e0);
    Real term1(0e0);
    Real term2(0e0);
    Real term(0e0);

    Real a_alpha_n(0e0);
    Real alpha_n(0e0);
    Real J0_r0_alpha_n(0e0);
    Real J1_a_alpha_n(0e0);

    Integer n(1);
    for(; n < num_term_use; ++n)
    {
        a_alpha_n = gsl_sf_bessel_zero_J0(n);
        alpha_n = a_alpha_n / a;
        J0_r0_alpha_n = gsl_sf_bessel_J0(r_0 * alpha_n);
        J1_a_alpha_n  = gsl_sf_bessel_J1(a_alpha_n);

        term1 = std::exp(-1e0 * alpha_n * alpha_n * Dt) * J0_r0_alpha_n;
        term2 = alpha_n * J1_a_alpha_n;
        term = term1 / term2;
        sum += term;

//             std::cout << "sum " << sum << ", term" << term << std::endl;

        if(fabs(term/sum) < threshold)
        {
//                 std::cout << "normal exit. " << n << std::endl;
            break;
        }
    }
    if(n == num_term_use)
        std::cout << "warning: use term over num_term_use" << std::endl;
    return (2e0 * sum / a);
}
Example #18
0
FELAnalysis::FELAnalysis(undulator &unduP, electronBeam &elecP, FELradiation &radiP)
{
    gamma0  = elecP.get_centralEnergy();
    sigmag0 = elecP.get_energySpread ();
    emitn   = elecP.get_emitnx       ();
    avgbeta = elecP.get_avgBetaFunc  ();
    lambdau = unduP.get_period       ();
    bfield  = unduP.get_field        ();
    current = elecP.get_peakCurrent  ();
    lambdas = radiP.get_wavelength   ();

    double b, JJ, etad, etae, etag, CapG;
    sigmax = sqrt(avgbeta*emitn/gamma0);
    au = sqrt(lambdas*2*gamma0*gamma0/lambdau-1);
    //au     = 0.934*bfield*lambdau*100/sqrt(2.0);
    b      = au*au/2.0/(1+au*au);
    JJ     = gsl_sf_bessel_J0(b) - gsl_sf_bessel_J1(b);
    rho1D  = pow((0.5/gamma0)*(0.5/gamma0)*(0.5/gamma0)*current/IA*(au*lambdau*JJ/2.0/PI/sigmax)*(au*lambdau*JJ/2.0/PI/sigmax),1.0/3.0);
    Lg1D   = lambdau/4/PI/sqrt(3)/rho1D;
    etad = Lg1D/(4*PI*sigmax*sigmax/lambdas);
    etae = Lg1D/avgbeta*4*PI*emitn/gamma0/lambdas;
    etag = Lg1D/lambdau*4*PI*sigmag0/gamma0;
    CapG = coefs.a1 *pow(etad,coefs.a2 )
         + coefs.a3 *pow(etae,coefs.a4 )
         + coefs.a5 *pow(etag,coefs.a6 )
         + coefs.a7 *pow(etae,coefs.a8 )*pow(etag,coefs.a9 )
         + coefs.a10*pow(etad,coefs.a11)*pow(etag,coefs.a12)
         + coefs.a13*pow(etad,coefs.a14)*pow(etae,coefs.a15)
         + coefs.a16*pow(etad,coefs.a17)*pow(etae,coefs.a18)*pow(etag,coefs.a19);
    Lg3D  = Lg1D*(1+CapG);
    rho3D = lambdau/4/PI/sqrt(3)/Lg3D;
    Psat  = 1.6*rho1D*Lg1D*Lg1D/(Lg3D*Lg3D)*gamma0*0.511*current/1000.0; // GW

    std::cout << "Lg1D : " << Lg1D  << " m"  << std::endl;
    std::cout << "Lg3D : " << Lg3D  << " m"  << std::endl;
    std::cout << "rho1D: " << rho1D          << std::endl;
    std::cout << "rho3D: " << rho3D          << std::endl;
    std::cout << "Psat : " << Psat  << " GW" << std::endl;
}
Example #19
0
void lamb_dipole(Vect3 centre, Array <Vect3> &X, Array <Vect3> &Omega, REAL amplitude, REAL theta_dipole, REAL radius, REAL a, REAL k) {
    //  Produce a Lamb Dipole, eg
    //        lamb_dipole(Vect3(100, 50, 0), X, OMEGA, 30, pi / 2, 30, 30, 3.81 / 30);


    REAL U;
    for (int x = 1; x <= NX; ++x) {
        for (int y = 1; y <= NY; ++y) {
            U = 0.;
            if (x < centre.x) U = -amplitude;
            if (x > centre.x) U = amplitude;

            REAL r = sqrt((centre[0] - x)*(centre[0] - x) + (centre[1] - y)*(centre[1] - y));
            if (r < radius) {
                REAL theta = atan2((centre[1] - y)*(centre[1] - y) + 1e-16, (centre[0] - x)*(centre[0] - x) + 1e-16);
                REAL numerator = gsl_sf_bessel_J1(k * r);
                REAL denominator = gsl_sf_bessel_J0(k * a);
                Omega.push_back(Vect3(0.0, 0.0, -U * 2 * k * sin(theta - theta_dipole) * numerator / denominator));
                X.push_back(Vect3((REAL) x, (REAL) y, 0.0));
            }
        }
    }
}
Example #20
0
static double bessel_J0(double w, void *p)
{
  return gsl_sf_bessel_J0(w);
}
Example #21
0
int main(int argc, char* argv[]) {
    printf("J0(%g) = %.18e\n", atof(argv[1]), gsl_sf_bessel_J0(atof(argv[1])));
    return 0;
}
Example #22
0
double F_born_wolf_psf_imag(double const rho, void *params)
{
    struct F_born_wolf_psf_params *p
        = (struct F_born_wolf_psf_params *) params;
    return gsl_sf_bessel_J0(p->alpha_r * rho) * -sin(p->psi * rho * rho) * rho;
}
Example #23
0
double FC_FUNC_(oct_bessel_j0, OCT_BESSEL_J0)
     (const double *x)
{
  return gsl_sf_bessel_J0(*x);
}
Example #24
0
void *writemod(void *work) {

SHARED_DATA *mydata = (SHARED_DATA *)work;


int nu0, nu1  ;
if (nui==-1){nu0=0;nu1=nnu[cIF];}else{nu0=nui;nu1=nui+1;};

int k,i,t,j,m,p, currow[6];

double phase, uu, vv, UVRad, Ampli, DerivRe[npar], DerivIm[npar];
double SPA, CPA, tA, tB, tempChi, deltat;
double wgtcorrA, tempRe[npar+1], tempIm[npar+1];
tempChi = 0.0;

const double deg2rad = 0.017453292519943295;

int Iam = mydata->Iam; 
int widx = (Iam == -1)?0:Iam;
int pmax = (mode == -1)?npar+1:1;
int mmax = (mode == 0)?1:ncomp;
bool write2mod = (mode == -3 || mode >= 0);
bool allmod = (mode == -3);
bool writeDer = (mode==-1);

double tempD, cosphase, sinphase, PA;
//double auxPhase, auxUVRad, auxPA;




for (t = mydata->t0[cIF]; t < mydata->t1[cIF] ; t++) {


 deltat = dt[cIF][t];

 if (fittable[cIF][t]!=0) {

  for (i = nu0; i < nu1; i++) {

    j = (nui!=-1)?0:i;
    k = nnu[cIF]*t+i;
    tempD = pow(wgt[0][cIF][k],2.);
 //   tempD = wgt[0][cIF][k];

    uu = freqs[cIF][i]*uv[0][cIF][t];
    vv =  freqs[cIF][i]*uv[1][cIF][t];

    for(p=0;p<pmax;p++){
      tempRe[p]=0.0; tempIm[p]=0.0;
    };

    for (m=0;m<mmax;m++){

   /*   auxUVRad = -1.0;
      auxPhase = 0.0;
      cosphase = 1.0;
      sinphase = 0.0;
      phase = 0.0;
      UVRad = 0.0;
      Ampli = 1.0;
      auxPA = 0.0;
      PA= 0.0; 
      SPA = 0.0; CPA = 1.0;  */

      for(p=0;p<6;p++){currow[p] = m*maxnchan*6+j+p*maxnchan;};
        wgtcorrA = exp(wgt[1][cIF][m*nt[cIF]+t]*pow(freqs[cIF][i],2.));

      for(p=0;p<pmax;p++){

        phase = (vars[p][currow[0]]+muRA[m]*deltat)*uu + (vars[p][currow[1]]+muDec[m]*deltat)*vv;
    //    if (phase != auxPhase){
    //      auxPhase = phase;
          cosphase = cos(phase);
          sinphase = sin(phase);
    //    };


        if (models[m] != 0) {
         PA = vars[p][currow[5]]*deg2rad;
    //     if (auxPA != PA){
           SPA = sin(PA);
           CPA = cos(PA);
     //      auxPA = PA;
     //    };
         tA = (uu*CPA - vv*SPA)*vars[p][currow[4]];
         tB = (uu*SPA + vv*CPA);
         UVRad = pow(tA*tA+tB*tB,0.5)*(vars[p][currow[3]]/2.0);



    //     if (UVRad != auxUVRad) {
    //      auxUVRad = UVRad;

         if (UVRad > 0.0) {
          switch (models[m]) {
            case 1: Ampli = (vars[p][currow[2]])*exp(-0.3606737602*UVRad*UVRad); break;
            case 2: Ampli = 2.0*(vars[p][currow[2]])*gsl_sf_bessel_J1(UVRad)/UVRad; break;
            case 3: Ampli = (vars[p][currow[2]])*gsl_sf_bessel_J0(UVRad); break;
            case 4: Ampli = 3.0*(vars[p][currow[2]])*(sin(UVRad)-UVRad*cos(UVRad))*(pow(UVRad,-3.0)); break;
            case 5: Ampli = (vars[p][currow[2]])*sin(UVRad)/UVRad; break;
            case 6: Ampli = (vars[p][currow[2]])*pow(1.+0.52034224525*UVRad*UVRad,-1.5); break;
            case 7: Ampli = 0.459224094*(vars[p][currow[2]])*gsl_sf_bessel_K0(UVRad); break;
            case 8: Ampli = (vars[p][currow[2]])*exp(-UVRad*1.3047660265); break;
            default: Ampli = vars[p][currow[2]];
          };
        } else {wgt[0][cIF][k]=0.0; Ampli=vars[p][currow[2]];};


    //    };

       } else {Ampli = vars[p][currow[2]];};

    Ampli *= wgtcorrA;

    tempRe[p] += Ampli*cosphase;
    tempIm[p] += Ampli*sinphase;

    if(allmod && m==0){
       ModVis[0][cIF][k] *= fixp[p][j];
       ModVis[1][cIF][k] *= fixp[p][j];
    };

    if(allmod || m==0) {
      if(write2mod){
         ModVis[0][cIF][k] += tempRe[p];
         ModVis[1][cIF][k] += tempIm[p];
      } else if(m==0) {
         tempRe[p] += ModVis[0][cIF][k]*fixp[p][j];
         tempIm[p] += ModVis[1][cIF][k]*fixp[p][j];
      };
    };

  };

  };

  tempChi += pow((ObsVis[0][cIF][k]-tempRe[0])*wgt[0][cIF][k],2.);
  tempChi += pow((ObsVis[1][cIF][k]-tempIm[0])*wgt[0][cIF][k],2.);


   if (writeDer) {

    for(p=0;p<npar;p++){
      DerivRe[p] = (tempRe[p+1]-tempRe[0])/dpar[p]; 
      DerivIm[p] = (tempIm[p+1]-tempIm[0])/dpar[p]; 

      WorkGrad[widx][p] += (ObsVis[0][cIF][k]-tempRe[0])*tempD*DerivRe[p];
      WorkGrad[widx][p] += (ObsVis[1][cIF][k]-tempIm[0])*tempD*DerivIm[p];
    };  

    for(p=0;p<npar;p++){
      for(m=p;m<npar;m++){
        WorkHess[widx][npar*p+m] += tempD*(DerivRe[p]*DerivRe[m] + DerivIm[p]*DerivIm[m]);
      };
    };

   };



  };

  };


};

if (writeDer){
for(p=0;p<npar;p++){
  for(m=p;m<npar;m++){
    WorkHess[widx][npar*m+p] = WorkHess[widx][npar*p+m];
  };
};
};

if (Iam != -1){Chi2[Iam] = tempChi; pthread_exit((void*) 0);} 
else {Chi2[0] = tempChi; return (void*) 0;};

}