int cosmic_ray_diffusion_evaluate(int target, int mode,
				  double *cr_E0_in, double *cr_E0_out, double *cr_E0_sum,
				  double *cr_n0_in, double *cr_n0_out, double *cr_n0_sum,
				  int *nexport, int *nsend_local, int CRpop)
{
  int startnode, numngb, listindex = 0;
  int j, n;
  MyDouble *pos;
  MyFloat h_i, rho;
  double dx, dy, dz;
  double h_i2, hinv_i, hinv4_i, hinv_j, hinv4_j;
  double dwk_i, h_j, dwk_j, dwk;
  double r, r2, u, CR_E0_Kappa_i, kappa_mean, w, wfac, cr_E0_out_sum, cr_E0_w_sum;
  double CR_n0_Kappa_i, cr_n0_out_sum, cr_n0_w_sum;

  if(mode == 0)
    {
      pos = P[target].Pos;
      h_i = PPP[target].Hsml;
      rho = SphP[target].d.Density;
      CR_E0_Kappa_i = CR_E0_Kappa[CRpop][target];
      CR_n0_Kappa_i = CR_n0_Kappa[CRpop][target];
    }
  else
    {
      pos = CR_DiffusionDataGet[target].Pos;
      h_i = CR_DiffusionDataGet[target].Hsml;
      rho = CR_DiffusionDataGet[target].Density;
      CR_E0_Kappa_i = CR_DiffusionDataGet[target].CR_E0_Kappa[CRpop];
      CR_n0_Kappa_i = CR_DiffusionDataGet[target].CR_n0_Kappa[CRpop];
    }

  h_i2 = h_i * h_i;
  hinv_i = 1.0 / h_i;
#ifndef  TWODIMS
  hinv4_i = hinv_i * hinv_i * hinv_i * hinv_i;
#else
  hinv4_i = hinv_i * hinv_i * hinv_i / boxSize_Z;
#endif

  /* initialize variables before SPH loop is started */
  cr_E0_out_sum = 0;
  cr_E0_w_sum = 0;
  cr_n0_out_sum = 0;
  cr_n0_w_sum = 0;

  /* Now start the actual SPH computation for this particle */

  if(mode == 0)
    {
      startnode = All.MaxPart;	/* root node */
    }
  else
    {
      startnode = CR_DiffusionDataGet[target].NodeList[0];
      startnode = Nodes[startnode].u.d.nextnode;	/* open it */
    }

  while(startnode >= 0)
    {
      while(startnode >= 0)
	{
	  numngb = ngb_treefind_pairs(pos, h_i, target, &startnode, mode, nexport, nsend_local);

	  if(numngb < 0)
	    return -1;

	  for(n = 0; n < numngb; n++)
	    {
	      j = Ngblist[n];

	      dx = pos[0] - P[j].Pos[0];
	      dy = pos[1] - P[j].Pos[1];
	      dz = pos[2] - P[j].Pos[2];

#ifdef PERIODIC			/*  now find the closest image in the given box size  */
	      if(dx > boxHalf_X)
		dx -= boxSize_X;
	      if(dx < -boxHalf_X)
		dx += boxSize_X;
	      if(dy > boxHalf_Y)
		dy -= boxSize_Y;
	      if(dy < -boxHalf_Y)
		dy += boxSize_Y;
	      if(dz > boxHalf_Z)
		dz -= boxSize_Z;
	      if(dz < -boxHalf_Z)
		dz += boxSize_Z;
#endif
	      r2 = dx * dx + dy * dy + dz * dz;

	      h_j = PPP[j].Hsml;
	      if(r2 < h_i2 || r2 < h_j * h_j)
		{
		  r = sqrt(r2);
		  if(r > 0)
		    {
		      if(r2 < h_i2)
			{
			  u = r * hinv_i;
			  if(u < 0.5)
			    dwk_i = hinv4_i * u * (KERNEL_COEFF_3 * u - KERNEL_COEFF_4);
			  else
			    dwk_i = hinv4_i * KERNEL_COEFF_6 * (1.0 - u) * (1.0 - u);
			}
		      else
			dwk_i = 0;

		      if(r2 < h_j * h_j)
			{
			  hinv_j = 1.0 / h_j;
#ifndef  TWODIMS
			  hinv4_j = hinv_j * hinv_j * hinv_j * hinv_j;
#else
			  hinv4_j = hinv_j * hinv_j * hinv_j / boxSize_Z;
#endif
			  u = r * hinv_j;
			  if(u < 0.5)
			    dwk_j = hinv4_j * u * (KERNEL_COEFF_3 * u - KERNEL_COEFF_4);
			  else
			    dwk_j = hinv4_j * KERNEL_COEFF_6 * (1.0 - u) * (1.0 - u);
			}
		      else
			dwk_j = 0;


		      dwk = 0.5 * (dwk_i + dwk_j);
		      wfac = 2.0 * P[j].Mass / (rho * SphP[j].d.Density) * (-dwk) / r;

		      /* cosmic ray diffusion equation kernel */
		      if((CR_E0_Kappa_i + CR_E0_Kappa[CRpop][j]) > 0)
			kappa_mean =
			  2 * (CR_E0_Kappa_i * CR_E0_Kappa[CRpop][j]) / (CR_E0_Kappa_i +
									 CR_E0_Kappa[CRpop][j]);
		      else
			kappa_mean = 0;

		      w = wfac * kappa_mean;
		      cr_E0_out_sum += (-w * cr_E0_in[j]);
		      cr_E0_w_sum += w;


		      /* cosmic ray diffusion equation kernel */
		      if((CR_n0_Kappa_i + CR_n0_Kappa[CRpop][j]) > 0)
			kappa_mean =
			  2 * (CR_n0_Kappa_i * CR_n0_Kappa[CRpop][j]) / (CR_n0_Kappa_i +
									 CR_n0_Kappa[CRpop][j]);
		      else
			kappa_mean = 0;

		      w = wfac * kappa_mean;
		      cr_n0_out_sum += (-w * cr_n0_in[j]);
		      cr_n0_w_sum += w;
		    }
		}
	    }
	}

      if(mode == 1)
	{
	  listindex++;
	  if(listindex < NODELISTLENGTH)
	    {
	      startnode = CR_DiffusionDataGet[target].NodeList[listindex];
	      if(startnode >= 0)
		startnode = Nodes[startnode].u.d.nextnode;	/* open it */
	    }
	}
    }


  /* Now collect the result at the right place */
  if(mode == 0)
    {
      cr_E0_out[target] = cr_E0_out_sum;
      cr_E0_sum[target] = cr_E0_w_sum;
      cr_n0_out[target] = cr_n0_out_sum;
      cr_n0_sum[target] = cr_n0_w_sum;
    }
  else
    {
      CR_DiffusionDataResult[target].CR_E0_Out = cr_E0_out_sum;
      CR_DiffusionDataResult[target].CR_E0_Sum = cr_E0_w_sum;
      CR_DiffusionDataResult[target].CR_n0_Out = cr_n0_out_sum;
      CR_DiffusionDataResult[target].CR_n0_Sum = cr_n0_w_sum;
    }

  return 0;
}
Esempio n. 2
0
File: hydra.c Progetto: Ingwar/amuse
/*! This function is the 'core' of the SPH force computation. A target
 *  particle is specified which may either be local, or reside in the
 *  communication buffer.
 */
void hydro_evaluate(int target, int mode)
{
  int j, k, n, timestep, startnode, numngb;
  FLOAT *pos, *vel;
  FLOAT mass, h_i, dhsmlDensityFactor, rho, pressure, f1, f2;
#ifdef MORRIS97VISC
  FLOAT alpha_visc, alpha_visc_j;
#endif
  double acc[3], dtEntropy, maxSignalVel;
  double dx, dy, dz, dvx, dvy, dvz;
  double h_i2, hinv, hinv4;
  double p_over_rho2_i, p_over_rho2_j, soundspeed_i, soundspeed_j;
#ifdef MONAGHAN83VISC
  double soundspeed_ij, h_ij;
#endif
  double hfc, dwk_i, vdotr, vdotr2, visc, mu_ij, rho_ij, vsig;
  double h_j, dwk_j, r, r2, u, hfc_visc;

#ifndef NOVISCOSITYLIMITER
  double dt;
#endif

  if(mode == 0)
    {
      pos = P[target].Pos;
      vel = SphP[target].VelPred;
      h_i = SphP[target].Hsml;
      mass = P[target].Mass;
      dhsmlDensityFactor = SphP[target].DhsmlDensityFactor;
      rho = SphP[target].Density;
      pressure = SphP[target].Pressure;
      timestep = P[target].Ti_endstep - P[target].Ti_begstep;
      soundspeed_i = sqrt(GAMMA * pressure / rho);
#ifdef MORRIS97VISC
      alpha_visc = SphP[target].Alpha;
#else
      f1 = fabs(SphP[target].DivVel) /
	(fabs(SphP[target].DivVel) + SphP[target].CurlVel +
	 0.0001 * soundspeed_i / SphP[target].Hsml / fac_mu);
#endif
    }
  else
    {
      pos = HydroDataGet[target].Pos;
      vel = HydroDataGet[target].Vel;
      h_i = HydroDataGet[target].Hsml;
      mass = HydroDataGet[target].Mass;
      dhsmlDensityFactor = HydroDataGet[target].DhsmlDensityFactor;
      rho = HydroDataGet[target].Density;
      pressure = HydroDataGet[target].Pressure;
      timestep = HydroDataGet[target].Timestep;
      soundspeed_i = sqrt(GAMMA * pressure / rho);
      f1 = HydroDataGet[target].F1;
#ifdef MORRIS97VISC
      alpha_visc = HydroDataGet[target].Alpha;
#endif
    }


  /* initialize variables before SPH loop is started */
  acc[0] = acc[1] = acc[2] = dtEntropy = 0;
  maxSignalVel = 0;

  p_over_rho2_i = pressure / (rho * rho) * dhsmlDensityFactor;
  h_i2 = h_i * h_i;

  /* Now start the actual SPH computation for this particle */
  startnode = All.MaxPart;
  do
    {
      numngb = ngb_treefind_pairs(&pos[0], h_i, &startnode);

      for(n = 0; n < numngb; n++)
	{
	  j = Ngblist[n];

	  dx = pos[0] - P[j].Pos[0];
	  dy = pos[1] - P[j].Pos[1];
	  dz = pos[2] - P[j].Pos[2];
#ifdef MORRIS97VISC
          alpha_visc_j = SphP[j].Alpha;
#endif	


#ifdef PERIODIC			/*  find the closest image in the given box size  */
	  if(dx > boxHalf_X)
	    dx -= boxSize_X;
	  if(dx < -boxHalf_X)
	    dx += boxSize_X;
	  if(dy > boxHalf_Y)
	    dy -= boxSize_Y;
	  if(dy < -boxHalf_Y)
	    dy += boxSize_Y;
	  if(dz > boxHalf_Z)
	    dz -= boxSize_Z;
	  if(dz < -boxHalf_Z)
	    dz += boxSize_Z;
#endif
	  r2 = dx * dx + dy * dy + dz * dz;
	  h_j = SphP[j].Hsml;
	  if(r2 < h_i2 || r2 < h_j * h_j)
	    {
	      r = sqrt(r2);
	      if(r > 0)
		{
		  p_over_rho2_j = SphP[j].Pressure / (SphP[j].Density * SphP[j].Density);
		  soundspeed_j = sqrt(GAMMA * p_over_rho2_j * SphP[j].Density);
		  dvx = vel[0] - SphP[j].VelPred[0];
		  dvy = vel[1] - SphP[j].VelPred[1];
		  dvz = vel[2] - SphP[j].VelPred[2];
		  vdotr = dx * dvx + dy * dvy + dz * dvz;

		  if(All.ComovingIntegrationOn)
		    vdotr2 = vdotr + hubble_a2 * r2;
		  else
		    vdotr2 = vdotr;

		  if(r2 < h_i2)
		    {
		      hinv = 1.0 / h_i;
#ifndef  TWODIMS
		      hinv4 = hinv * hinv * hinv * hinv;
#else
		      hinv4 = hinv * hinv * hinv / boxSize_Z;
#endif
		      u = r * hinv;
		      if(u < 0.5)
			dwk_i = hinv4 * u * (KERNEL_COEFF_3 * u - KERNEL_COEFF_4);
		      else
			dwk_i = hinv4 * KERNEL_COEFF_6 * (1.0 - u) * (1.0 - u);
		    }
		  else
		    {
		      dwk_i = 0;
		    }

		  if(r2 < h_j * h_j)
		    {
		      hinv = 1.0 / h_j;
#ifndef  TWODIMS
		      hinv4 = hinv * hinv * hinv * hinv;
#else
		      hinv4 = hinv * hinv * hinv / boxSize_Z;
#endif
		      u = r * hinv;
		      if(u < 0.5)
			dwk_j = hinv4 * u * (KERNEL_COEFF_3 * u - KERNEL_COEFF_4);
		      else
			dwk_j = hinv4 * KERNEL_COEFF_6 * (1.0 - u) * (1.0 - u);
		    }
		  else
		    {
		      dwk_j = 0;
		    }

		  if(soundspeed_i + soundspeed_j > maxSignalVel)
		    maxSignalVel = soundspeed_i + soundspeed_j;

		  if(vdotr2 < 0)	/* ... artificial viscosity */
		    {
#ifndef MONAGHAN83VISC
		      mu_ij = fac_mu * vdotr2 / r;	/* note: this is negative! */
#else		      
		      h_ij = 0.5 * (h_i + h_j);
		      mu_ij = fac_mu * h_ij * vdotr2 / (r2 + 0.0001 * h_ij * h_ij);
#endif
                      vsig = soundspeed_i + soundspeed_j - 3 * mu_ij;



		      if(vsig > maxSignalVel)
			maxSignalVel = vsig;

		      rho_ij = 0.5 * (rho + SphP[j].Density);

#ifdef MORRIS97VISC		      
                      visc = 0.25 * (alpha_visc + alpha_visc_j) * vsig * (-mu_ij) / rho_ij;
#else

		      f2 =
			fabs(SphP[j].DivVel) / (fabs(SphP[j].DivVel) + SphP[j].CurlVel +
						0.0001 * soundspeed_j / fac_mu / SphP[j].Hsml);

#ifndef MONAGHAN83VISC
		      visc = 0.25 * All.ArtBulkViscConst * vsig * (-mu_ij) / rho_ij * (f1 + f2);
#else			      
		      soundspeed_ij = (soundspeed_i+soundspeed_j) * 0.5;

		      visc = ((-All.ArtBulkViscConst) * soundspeed_ij * mu_ij + All.ArtBulkViscBeta * mu_ij * mu_ij) / rho_ij;
#endif //MONAGHAN83VISC

#endif //MORRIS97VISC

		      /* .... end artificial viscosity evaluation */
#ifndef NOVISCOSITYLIMITER
		      /* make sure that viscous acceleration is not too large */
		      dt = imax(timestep, (P[j].Ti_endstep - P[j].Ti_begstep)) * All.Timebase_interval;
		      if(dt > 0 && (dwk_i + dwk_j) < 0)
			{
			  visc = dmin(visc, 0.5 * fac_vsic_fix * vdotr2 /
				      (0.5 * (mass + P[j].Mass) * (dwk_i + dwk_j) * r * dt));
			}
#endif
		    }
		  else
		    visc = 0;

		  p_over_rho2_j *= SphP[j].DhsmlDensityFactor;

		  hfc_visc = 0.5 * P[j].Mass * visc * (dwk_i + dwk_j) / r;

		  hfc = hfc_visc + P[j].Mass * (p_over_rho2_i * dwk_i + p_over_rho2_j * dwk_j) / r;

		  acc[0] -= hfc * dx;
		  acc[1] -= hfc * dy;
		  acc[2] -= hfc * dz;
		  dtEntropy += 0.5 * hfc_visc * vdotr2;
		}
	    }
	}
    }
  while(startnode >= 0);

  /* Now collect the result at the right place */
  if(mode == 0)
    {
      for(k = 0; k < 3; k++)
	SphP[target].HydroAccel[k] = acc[k];
      SphP[target].DtEntropy = dtEntropy;
      SphP[target].MaxSignalVel = maxSignalVel;
    }
  else
    {
      for(k = 0; k < 3; k++)
	HydroDataResult[target].Acc[k] = acc[k];
      HydroDataResult[target].DtEntropy = dtEntropy;
      HydroDataResult[target].MaxSignalVel = maxSignalVel;
    }
}
Esempio n. 3
0
/* this function evaluates parts of the matrix A */
int radtransfer_evaluate(int target, int mode, double *in, double *out, double *sum, int *nexport,
			 int *nsend_local)
{
  int startnode, numngb, listindex = 0;
  int j, n, k;
  MyFloat *ET_aux, ET_j[6], ET_i[6], ET_ij[6];
  MyFloat kappa_i, kappa_j, kappa_ij;
#ifdef RADTRANSFER_FLUXLIMITER
  MyFloat lambda_i, lambda_j;
#endif
  MyDouble *pos;
  MyFloat mass, mass_i, rho, rho_i;
  double sum_out = 0, sum_w = 0, a3, fac = 0;

  double dx, dy, dz;
  double h_j, hinv, hinv4, h_i;
  double dwk_i, dwk_j, dwk;
  double r, r2, r3, u;

#ifdef PERIODIC
  double boxsize, boxhalf;

  boxsize = All.BoxSize;
  boxhalf = 0.5 * All.BoxSize;
#endif

  if(All.ComovingIntegrationOn)
    a3 = All.Time * All.Time * All.Time;
  else
    a3 = 1.0;

  if(mode == 0)
    {
      ET_aux = SphP[target].ET;
      pos = P[target].Pos;
      h_i = PPP[target].Hsml;
      kappa_i = Kappa[target];
#ifdef RADTRANSFER_FLUXLIMITER
      lambda_i = Lambda[target];
#endif
      mass_i = P[target].Mass;
      rho_i = SphP[target].d.Density;
    }
  else
    {
      ET_aux = RadTransferDataGet[target].ET;
      pos = RadTransferDataGet[target].Pos;
      h_i = RadTransferDataGet[target].Hsml;
      kappa_i = RadTransferDataGet[target].Kappa;
#ifdef RADTRANSFER_FLUXLIMITER
      lambda_i = RadTransferDataGet[target].Lambda;
#endif
      mass_i = RadTransferDataGet[target].Mass;
      rho_i = RadTransferDataGet[target].Density;
    }

#ifdef RADTRANSFER_MODIFY_EDDINGTON_TENSOR
  /*modify Eddington tensor */
  ET_i[0] = 2 * ET_aux[0] - 0.5 * ET_aux[1] - 0.5 * ET_aux[2];
  ET_i[1] = 2 * ET_aux[1] - 0.5 * ET_aux[2] - 0.5 * ET_aux[0];
  ET_i[2] = 2 * ET_aux[2] - 0.5 * ET_aux[0] - 0.5 * ET_aux[1];

  for(k = 3; k < 6; k++)
    ET_i[k] = 2.5 * ET_aux[k];
#else
  for(k = 0; k < 6; k++)
    ET_i[k] = ET_aux[k];
#endif

  if(mode == 0)
    {
      startnode = All.MaxPart;
    }
  else
    {
      startnode = RadTransferDataGet[target].NodeList[0];
      startnode = Nodes[startnode].u.d.nextnode;
    }

  while(startnode >= 0)
    {
      while(startnode >= 0)
	{
	  numngb = ngb_treefind_pairs(pos, h_i, target, &startnode, mode, nexport, nsend_local);

	  if(numngb < 0)
	    return -1;

	  for(n = 0; n < numngb; n++)
	    {
	      j = Ngblist[n];
#ifdef RT_VAR_DT
	      if(SphP[j].rt_flag == 1)
#endif
		{
		  dx = pos[0] - P[j].Pos[0];
		  dy = pos[1] - P[j].Pos[1];
		  dz = pos[2] - P[j].Pos[2];
#ifdef PERIODIC			/*  now find the closest image in the given box size  */
		  if(dx > boxHalf_X)
		    dx -= boxSize_X;
		  if(dx < -boxHalf_X)
		    dx += boxSize_X;
		  if(dy > boxHalf_Y)
		    dy -= boxSize_Y;
		  if(dy < -boxHalf_Y)
		    dy += boxSize_Y;
		  if(dz > boxHalf_Z)
		    dz -= boxSize_Z;
		  if(dz < -boxHalf_Z)
		    dz += boxSize_Z;
#endif
		  r2 = dx * dx + dy * dy + dz * dz;
		  r = sqrt(r2);
		  r3 = r2 * r;
		  h_j = PPP[j].Hsml;
		  
		  if(r > 0 && (r < h_i || r < h_j))
		    {
		      mass = P[j].Mass;
		      rho = SphP[j].d.Density;
		      kappa_j = Kappa[j];
#ifdef RADTRANSFER_FLUXLIMITER
		      lambda_j = Lambda[j];
#endif
		      
#ifdef RADTRANSFER_MODIFY_EDDINGTON_TENSOR
		      ET_aux = SphP[j].ET;
		      
		      /*modify Eddington tensor */
		      ET_j[0] = 2 * ET_aux[0] - 0.5 * ET_aux[1] - 0.5 * ET_aux[2];
		      ET_j[1] = 2 * ET_aux[1] - 0.5 * ET_aux[2] - 0.5 * ET_aux[0];
		      ET_j[2] = 2 * ET_aux[2] - 0.5 * ET_aux[0] - 0.5 * ET_aux[1];
		      
		      for(k = 3; k < 6; k++)
			ET_j[k] = 2.5 * ET_aux[k];
#else
		      for(k = 0; k < 6; k++)
			ET_j[k] = SphP[j].ET[k];
#endif
		      
		      for(k = 0; k < 6; k++)
			ET_ij[k] = 0.5 * (ET_i[k] + ET_j[k]);
		      
		      if(r < h_i)
			{
			  hinv = 1.0 / h_i;
			  hinv4 = hinv * hinv * hinv * hinv;
			  u = r * hinv;
			  
			  if(u < 0.5)
			    dwk_i = hinv4 * u * (KERNEL_COEFF_3 * u - KERNEL_COEFF_4);
			  else
			    dwk_i = hinv4 * KERNEL_COEFF_6 * (1.0 - u) * (1.0 - u);
			}
		      else
			dwk_i = 0;
		      
		      if(r < h_j)
			{
			  hinv = 1.0 / h_j;
			  hinv4 = hinv * hinv * hinv * hinv;
			  u = r * hinv;
			  
			  if(u < 0.5)
			    dwk_j = hinv4 * u * (KERNEL_COEFF_3 * u - KERNEL_COEFF_4);
			  else
			    dwk_j = hinv4 * KERNEL_COEFF_6 * (1.0 - u) * (1.0 - u);
			}
		      else
			dwk_j = 0;
		      
		      kappa_ij = 0.5 * (1/kappa_i + 1/kappa_j);
		      dwk = 0.5 * (dwk_i + dwk_j);	 
		      mass = 0.5 * (mass + mass_i);
		      rho = 0.5 * (rho + rho_i);
		      
		      double tensor = (ET_ij[0] * dx * dx + ET_ij[1] * dy * dy + ET_ij[2] * dz * dz
				       + 2.0 * ET_ij[3] * dx * dy + 2.0 * ET_ij[4] * dy * dz + 2.0 * ET_ij[5] * dz * dx);
		      
		      if(tensor > 0)
			{
			  /* divide c_light by a because kappa is comoving */
			  /* all variables are comoving!!! */
			  fac = -2.0 * dt * (c_light / a3) * (mass / rho) * kappa_ij * dwk / r3 * tensor;
			  
#ifdef RADTRANSFER_FLUXLIMITER
			  fac *= 0.5 * (lambda_i + lambda_j);
#endif
			  
			  sum_out -= fac * in[j];
			  
			  sum_w += fac;
			}		
		    }
		}
	    }
	}
      
      if(mode == 1)
	{
	  listindex++;
	  if(listindex < NODELISTLENGTH)
	    {
	      startnode = RadTransferDataGet[target].NodeList[listindex];
	      if(startnode >= 0)
		startnode = Nodes[startnode].u.d.nextnode;
	    }
	}
    }
  
  if(mode == 0)
    {
      out[target] = sum_out;
      sum[target] = sum_w;
    }
  else
    {
      RadTransferDataResult[target].Out = sum_out;
      RadTransferDataResult[target].Sum = sum_w;
    }
  
  return 0;
}