Example #1
0
void run_output() {
	int i, j;
	int ipart;
	char filename[128];
	double pos[nDim];
	double phi;
	double x_rms, x_norm;
	double v_rms, v_norm;
	double x_rms_total, v_rms_total;
	double x_norm_total, v_norm_total;
	double x_z, v_z;
	double q, q_i;
	int num_level_cells;
	int *level_cells;

	sprintf(filename, "dumps/zeldovich_[%5.4f]_%03d.dat", abox[min_level], local_proc_id );
	output = fopen(filename, "w");

#ifdef PARTICLES
	sprintf(filename, "dumps/zeldovich_particles_[%5.4f]_%03d.dat", abox[min_level], local_proc_id );
	particles = fopen(filename, "w");
#endif
	dgrowth0 = growth(auni_init);
	dgrowth = growth(abox[min_level]);
        ddgrowthdt = dgrowthdt(abox[min_level]);

	c1 = 6.0/abox[min_level];
	c2 = 1.0 / (ak*ak);
	c3 = dgrowth*ampl*ak;
#ifdef GRAVITY
	pos[1] = pos[2] = (double)num_grid / 2.0 - 0.5*cell_size[max_level];
	pos[0] = 0.0;
	i = cell_find_position(pos);
	if ( i != -1 ) {
		phi = cell_potential(i);
	} else { 
		phi = 0.0;
	}
#endif

	MPI_Allreduce( &phi, &phi0, 1, MPI_DOUBLE, MPI_MAX, mpi.comm.run );

	/* uncommment to print all cells */
/*
	for ( i = min_level; i <= max_level; i++ ) {
		select_level( i, CELL_TYPE_LOCAL, &num_level_cells, &level_cells );
		for ( j = 0; j < num_level_cells; j++ ) {
			output_cell( level_cells[j], i );
		}
		cart_free(level_cells);
	}
*/

	for ( pos[0] = 0.0; pos[0] < num_grid; pos[0] += 0.5*cell_size[max_level] ) {
            i=cell_find_position(pos);
		if ( i != -1 ) {
			output_cell( i, cell_level(i) );
		}
	}

	fclose(output);

	/* compute rms fluctuations */
	x_rms = x_norm = 0.0;
	v_rms = v_norm = 0.0;

	ddgrowthdt = dgrowthdt( abox_from_tcode( tl[min_level] - 0.5*dtl[min_level] ) );

#ifdef PARTICLES
	for ( i = 0; i < num_particles; i++ ) {
		if ( particle_level[i] != FREE_PARTICLE_LEVEL ) {
			ipart = i;
			x0 = particle_x[i][0];
			q = root_finder( qsolve, 0.0, num_grid, 1e-9, 1e-9 );

			q_i = particle_q_init( particle_id[i] );

			x_z = q_i + dgrowth*ampl*sin(ak*q_i);
			v_z = ddgrowthdt*ampl*sin(ak*q_i);

			x_rms += (x0-x_z)*(x0-x_z);
			x_norm += (x_z-q_i)*(x_z-q_i);

			v_rms += (particle_v[i][0]-v_z)*(particle_v[i][0]-v_z);
			v_norm += v_z*v_z;
                        
			fprintf(particles, "%lu %e %e %e %e %e \n",
					particle_id[i],
					particle_x[i][0],
					particle_v[i][0], 
					particle_pot[i],
					x_z,v_z ); 
                        
                        
		}
	}
        
	fclose(particles);
	MPI_Reduce( &x_rms, &x_rms_total, 1, MPI_DOUBLE, MPI_SUM, MASTER_NODE, mpi.comm.run );
	MPI_Reduce( &v_rms, &v_rms_total, 1, MPI_DOUBLE, MPI_SUM, MASTER_NODE, mpi.comm.run );
	MPI_Reduce( &x_norm, &x_norm_total, 1, MPI_DOUBLE, MPI_SUM, MASTER_NODE, mpi.comm.run );
	MPI_Reduce( &v_norm, &v_norm_total, 1, MPI_DOUBLE, MPI_SUM, MASTER_NODE, mpi.comm.run );

	if ( local_proc_id == MASTER_NODE ) {	

		x_rms_total /= x_norm_total;
		v_rms_total /= v_norm_total;

		x_rms_total = sqrt(x_rms_total);
		v_rms_total = sqrt(v_rms_total);

		particles = fopen("dumps/particle_rms.dat", "a");
		fprintf(particles, "%e %e %e %e %e\n", abox[min_level], x_rms_total, v_rms_total, x_norm_total, v_norm_total );
		fclose(particles);
	}
#endif /* PARTICLES */
}
Example #2
0
void init_run() {
	int i, j, k;
	int index;
	double a_th;

	int ipart;
	int icell;
	double qi, qj, qk;
	double xcons, vcons;
	double dx, dvx;
	double pw;
	double a_vel;
	double qfact;

	int num_level_cells;
	int *level_cells;

	cosmology_set(OmegaM,OmM0);
	cosmology_set(OmegaB,OmB0);
	cosmology_set(OmegaL,OmL0);
	cosmology_set(h,h0);
	cosmology_set(DeltaDC,dDC);
	box_size = Lbox0;

	units_set_art(cosmology->OmegaM,cosmology->h,box_size);
	units_init();
	build_cell_buffer();
	repair_neighbors();
        
	auni[min_level] = auni_init;
	tl[min_level] = tcode_from_auni( auni_init );
	for ( i = min_level; i <= max_level; i++ ) { tl[i] = tl[min_level]; }
	abox[min_level] = auni_init;

	for(i=min_level+1; i<=max_level; i++)
	{
            tl[i] = tl[min_level];
            auni[i] = auni[min_level];
            abox[i] = abox[min_level];
	}
        
	units_update(min_level);
	cart_debug("tl[min_level] = %f", tl[min_level] );
	cart_debug("au[min_level] = %f", auni[min_level] );
	cart_debug("ab[min_level] = %f", abox[min_level] );
	cart_debug("DC mode = %f", cosmology->DeltaDC );
	cosmology_set_fixed();


	rhogas0 = cosmology->OmegaB/cosmology->OmegaM;
	cart_debug("rhogas0 = %e", rhogas0 );

	Tinit = TinitK/units->temperature;


	ak = 2.0*M_PI / lambda;
	dgrowth = growth(abox[min_level]);
	ddgrowthdt = dgrowthdt(abox[min_level]);
	ampl = 1.0 / ( growth(a_cross) * ak );
	cart_debug("Tinit,TinitK = %e %e", Tinit,TinitK );

#ifdef HYDRO
	for ( i = min_level; i <= max_level; i++ ) {
		cart_debug("generating initial conditions on level %u", i );
	
		select_level( i, CELL_TYPE_LOCAL, &num_level_cells, &level_cells );
		for ( j = 0; j < num_level_cells; j++ ) {
//                    cart_debug("%d %d",level_cells[j],num_cells);
			initial_conditions( level_cells[j], i );
		}
		cart_free( level_cells );
	}

	for ( i = min_level; i <= max_level; i++ ) {
		update_buffer_level( i, all_hydro_vars, num_hydro_vars );
	}
#endif /* HYDRO */

	cart_debug("choose timestep and set velocity on the half step");
	dtl[min_level] = 0.0;
	set_timestepping_scheme();
	dtl[min_level]=.125;
	cart_debug("=======================%e",dtl[min_level]);

	dtl_old[min_level] = dtl[min_level];
	tl_old[min_level] = tl[min_level]-dtl[min_level];
	abox_old[min_level] = abox_from_tcode(tl_old[min_level]);
	dtl_old[min_level] = dtl[min_level];

	for ( i = min_level+1; i <= max_level; i++ ) {
		tl_old[i] = tl[i]-dtl[i];
		abox_old[i] = abox_from_tcode(tl_old[i]);
		dtl_old[i] = dtl[i];
	}

#ifdef GRAVITY
#ifdef HYDRO
	for ( i = min_level; i <= max_level; i++ ) {
		cart_debug("generating gravity on level %u", i );
	
//                cart_assert(dtl[i]==0);
		select_level( i, CELL_TYPE_LOCAL, &num_level_cells, &level_cells );
		for ( j = 0; j < num_level_cells; j++ ) {
			initial_gravity( level_cells[j], i );
		}
		cart_free( level_cells );
	}

	for ( i = min_level; i <= max_level; i++ ) {
		update_buffer_level( i, all_hydro_vars, num_hydro_vars );
	}
#endif /* GRAVITY */
#endif /* HYDRO */
        
        
#ifdef PARTICLES
	qfact = (double)num_grid / (double)num_grid;
	pw = (1.0-rhogas0)*qfact*qfact*qfact;

	cart_debug("particle weight = %e", pw );

	xcons = dgrowth*ampl;
	a_vel = abox_from_tcode( tl[min_level] - 0.5*dtl[min_level] );
	vcons = ampl * dgrowthdt( a_vel);

	ipart = 0;
	for ( i = 0; i < num_grid; i++ ) {
		qi = qfact*((double)i + 0.5);
		dx = xcons * sin( ak * qi );
		dvx = vcons * sin( ak * qi );
		for ( j = 0; j < num_grid; j++ ) {
			qj = qfact*((double)j + 0.5);
			for ( k = 0; k < num_grid; k++ ) {
				qk = qfact*((double)k + 0.5);

				particle_x[ipart][0] = qi + dx;
				particle_x[ipart][1] = qj;
				particle_x[ipart][2] = qk;

				if ( particle_x[ipart][0] >= (double)num_grid ) {
					particle_x[ipart][0] -= num_grid;
				}

				if ( particle_x[ipart][1] >= (double)num_grid ) {
					particle_x[ipart][1] -= num_grid;
				}

				if ( particle_x[ipart][2] >= (double)num_grid ) {
					particle_x[ipart][2] -= num_grid;
				}

				icell = cell_find_position( particle_x[ipart] );

				if ( icell != -1 && cell_is_local(icell) ) {
					particle_v[ipart][0] = dvx;
					particle_v[ipart][1] = 0.0;
					particle_v[ipart][2] = 0.0;

					particle_id[ipart] = (particleid_t)num_grid*(num_grid*i + j) + k;
					particle_mass[ipart] = pw;

					cart_assert( qi == particle_q_init( particle_id[ipart] ) );

					particle_t[ipart] = tl[min_level];
					particle_dt[ipart] = dtl[min_level];

					ipart++;
				}
			}
		}
	}

	cart_debug("created %u particles", ipart );

	num_local_particles = ipart;
	num_particles_total = (particleid_t)num_grid*(particleid_t)num_grid*(particleid_t)num_grid;
	num_particle_species = 1;
	particle_species_mass[0] = pw;
	particle_species_num[0] = num_particles_total;
	particle_species_indices[0] = 0;
	particle_species_indices[1] = num_particles_total;

	build_particle_list();

/* 	assign_density( min_level, min_level ); */ //for refinement
/* 	modify( min_level, 0 ); */
 	assign_density( min_level, min_level );  //for refinement
 	modify( min_level, 0 ); 

	if ( local_proc_id == MASTER_NODE ) {
		particles = fopen("dumps/particle_rms.dat", "w");
		fclose(particles);
	}
#endif

	check_map();
	cart_debug("done with initialization");
}
Example #3
0
int main(int argc, char * argv[]) {


  fftw_plan pc2r;
  fftw_plan pr2c;
  fftw_complex *dvdk;
  DIR* dir;
  float *temp;
  char fname[256];
  FILE *fid,*fidtxt, *fidtxt2;
  double *t21, *dvdr; 
  long int i,j,p,indi,indj;
  int nmaxcut,nmincut;
  double aver,Tcmb,kk,growth,dgdt,Hubz,TS,xtot;
  double zmin,zmax,dz,z;


  /* Check for correct number of parameters*/
  if(argc != 2) {
    printf("Usage : t21 base_dir\n");
    exit(1);
  }
  get_Simfast21_params(argv[1]);
  if(global_use_Lya_xrays==0) printf("Lya and xray use set to false - assuming TS>>TCMB in t21 calculation.\n");
  zmin=global_Zminsim;
  zmax=global_Zmaxsim;
  dz=global_Dzsim;

#ifdef _OMPTHREAD_
  omp_set_num_threads(global_nthreads);
  fftw_init_threads();
  fftw_plan_with_nthreads(global_nthreads);
  printf("Using %d threads\n",global_nthreads);
#endif
 

  /* Create directory deltaTb */
  sprintf(fname,"%s/deltaTb",argv[1]);
  if((dir=opendir(fname))==NULL) {  
    printf("Creating t21 directory\n");
    if(mkdir(fname,(S_IRWXU | S_IRGRP | S_IXGRP | S_IROTH | S_IXOTH))!=0) {
      printf("Error creating directory!\n");
      exit(1);
    }
  }

  sprintf(fname,"%s/Output_text_files/t21_av_N%ld_L%.1f.dat",argv[1],global_N_smooth,global_L/global_hubble); 
  if((fidtxt=fopen(fname,"a"))==NULL){
    printf("\nError opening output %s file...\n",fname); 
    exit(1);
  }  
  sprintf(fname,"%s/Output_text_files/TS_av_N%ld_L%.1f.dat",argv[1],global_N_smooth,global_L/global_hubble); 
  if((fidtxt2=fopen(fname,"a"))==NULL){
    printf("\nError opening output %s file...\n",fname); 
    exit(1);
  }  
  if(!(temp=(float *) malloc(global_N3_smooth*sizeof(float)))) {
    printf("Problem...\n");
    exit(1);
  }
  if(!(t21=(double *) malloc(global_N3_smooth*sizeof(double)))) {
    printf("Problem...\n");
    exit(1);
  }
  if(!(dvdk = (fftw_complex*) fftw_malloc(sizeof(fftw_complex)*global_N_smooth*global_N_smooth*(global_N_smooth/2+1)))) {  
    printf("Problem...\n");
    exit(1);
  } 
  if(!(dvdr=(double *) fftw_malloc(global_N3_smooth*sizeof(double)))) {  
    printf("Problem...\n");
    exit(1);
  }
  if(!(pr2c=fftw_plan_dft_r2c_3d(global_N_smooth, global_N_smooth, global_N_smooth, dvdr , dvdk, FFTW_MEASURE))) { 
    printf("Problem...\n");
    exit(1);
  }
  if(!(pc2r=fftw_plan_dft_c2r_3d(global_N_smooth, global_N_smooth, global_N_smooth, dvdk, dvdr, FFTW_MEASURE))) { 
    printf("Problem...\n");
    exit(1);
  }

  /**************************************************/
  /************** redshift cycle ********************/
  /**************************************************/
  for(z=zmax;z>(zmin-dz/10);z-=dz){    
    printf("z = %f\n",z);fflush(0);
    sprintf(fname,"%s/deltaTb/deltaTb_z%.3lf_N%ld_L%.1f.dat",argv[1],z,global_N_smooth,global_L/global_hubble);
    if((fid = fopen(fname,"rb"))!=NULL) {
      printf("File:%s already exists - skipping this redshift...\n",fname);
      fclose(fid);
    }else {      
      if(z>(global_Zminsfr-global_Dzsim/10) && global_use_Lya_xrays==1) {
	Tcmb=Tcmb0*(1.+z);
	sprintf(fname,"%s/x_c/xc_z%.3lf_N%ld_L%.1f.dat",argv[1],z,global_N_smooth,global_L/global_hubble);
	if((fid = fopen(fname,"rb"))==NULL) {
	  printf("Error opening file:%s\n",fname);
	  exit(1);
	}
	fread(temp,sizeof(float),global_N3_smooth,fid);
	fclose(fid);
	for(i=0;i<global_N3_smooth;i++) t21[i]=(double)temp[i];
	sprintf(fname,"%s/Lya/xalpha_z%.3lf_N%ld_L%.1f.dat",argv[1],z,global_N_smooth,global_L/global_hubble);
	if((fid = fopen(fname,"rb"))==NULL) {
	  printf("Error opening file:%s\n",fname);
	  exit(1);
	}
	fread(temp,sizeof(float),global_N3_smooth,fid);
	fclose(fid);
	aver=0.;
	for(i=0;i<global_N3_smooth;i++) {
	  xtot=(double)temp[i]+t21[i];
	  t21[i]=xtot/(1.+xtot);
	}
	
	sprintf(fname,"%s/xrays/TempX_z%.3lf_N%ld_L%.1f.dat",argv[1],z,global_N_smooth,global_L/global_hubble);
	if((fid = fopen(fname,"rb"))==NULL) {
	  printf("Error opening file:%s\n",fname);
	  exit(1);
	}

	fread(temp,sizeof(float),global_N3_smooth,fid);
	fclose(fid);
	
	TS=0.;
	for(i=0;i<global_N3_smooth;i++) {
	
	  t21[i]=t21[i]*(1.-Tcmb/(double)temp[i]); // Temperature correction for high redshifts
	  TS+=Tcmb/(1.-t21[i]);
	}
      }else {
	for(i=0;i<global_N3_smooth;i++) t21[i]=1.0;
      }
      sprintf(fname, "%s/delta/deltanl_z%.3f_N%ld_L%.1f.dat",argv[1],z,global_N_smooth,global_L/global_hubble); 
      fid=fopen(fname,"rb");
      if (fid==NULL) {printf("Error reading deltanl file... Check path or if the file exists..."); exit (1);}
      fread(temp,sizeof(float),global_N3_smooth,fid);
      fclose(fid);
      for(i=0;i<global_N3_smooth;i++) dvdr[i]=(double)temp[i];
      /* Executes FFT */
      fftw_execute(pr2c);
      /********************************************************************/
      printf("Computing v field...\n"); 
#ifdef _OMPTHREAD_
#pragma omp parallel for shared(dvdk,global_dk) private(i,indi,j,indj,p,kk) 
#endif
      for(i=0;i<global_N_smooth;i++) {         
	if(i>global_N_smooth/2) {    /* Large frequencies are equivalent to smaller negative ones */
	  indi=-(global_N_smooth-i);
	}else indi=i;
	for(j=0;j<global_N_smooth;j++) {           
	  if(j>global_N_smooth/2) {  
	    indj=-(global_N_smooth-j);
	  }else indj=j;  
	  for(p=0;p<=global_N_smooth/2;p++) {
	    if(i==0 && j==0 && p==0) {
	      dvdk[i*global_N_smooth*(global_N_smooth/2+1)+j*(global_N_smooth/2+1)+p]=0.;
	    }else {
	      kk=global_dk*sqrt(indi*indi+indj*indj+p*p);	
	      if(global_vi==1) {
		dvdk[i*global_N_smooth*(global_N_smooth/2+1)+j*(global_N_smooth/2+1)+p]=(global_dk*indi)*(global_dk*indi)/(kk*kk)*dvdk[i*global_N_smooth*(global_N_smooth/2+1)+j*(global_N_smooth/2+1)+p];
	      }else if(global_vi==2) {
		dvdk[i*global_N_smooth*(global_N_smooth/2+1)+j*(global_N_smooth/2+1)+p]=(global_dk*indj)*(global_dk*indj)/(kk*kk)*dvdk[i*global_N_smooth*(global_N_smooth/2+1)+j*(global_N_smooth/2+1)+p];
	      }else {
		dvdk[i*global_N_smooth*(global_N_smooth/2+1)+j*(global_N_smooth/2+1)+p]=(global_dk*p)*(global_dk*p)/(kk*kk)*dvdk[i*global_N_smooth*(global_N_smooth/2+1)+j*(global_N_smooth/2+1)+p];
	      }
	    }
	  }
	}
      }
      /* Executes FFT */
      fftw_execute(pc2r);
      nmaxcut=0;
      nmincut=0;
      growth=getGrowth(z);
      dgdt=dgrowthdt(z);
      Hubz=Hz(z);
#ifdef _OMPTHREAD_
#pragma omp parallel for shared(dvdr,global_L,global_dx_smooth,global_N3_smooth) private(i) 
#endif
      for(i=0;i<global_N3_smooth;i++) {
	dvdr[i]=-dgdt*dvdr[i]/global_L/global_L/global_L/Hubz*global_dx_smooth*global_dx_smooth*global_dx_smooth/growth; /* FT normalization + divide by H(z). Divide by growth to get deltanl back to z=0 */
      }
      for(i=0;i<global_N3_smooth;i++) {
	if(dvdr[i] > maxcut) {dvdr[i]=maxcut; nmaxcut++;}
	else if(dvdr[i] < mincut) {dvdr[i]=mincut; nmincut++;}
      }
      printf("nmaxcut: %d (%f %%), nmincut: %d (%f %%)\n\n",nmaxcut,100.*nmaxcut/global_N3_smooth,nmincut,100.*nmincut/global_N3_smooth);fflush(0);
     
#ifdef _OMPTHREAD_
#pragma omp parallel for shared(t21,temp,dvdr,global_hubble,global_omega_b,global_omega_m,global_N3_smooth,z) private(i) 
#endif
      for(i=0;i<global_N3_smooth;i++) {
  	/* output in K */
	t21[i]=23.0/1000.*(1.+(double)temp[i])*t21[i]/(1.+1.*dvdr[i])*(0.7/global_hubble)*(global_omega_b*global_hubble*global_hubble/0.02)*sqrt((0.15/global_omega_m/global_hubble/global_hubble)*(1.+z)/10.);  /*Units in Kelvin */
      }
      
      sprintf(fname,"%s/Ionization/xHII_z%.3f_eff%.2lf_N%ld_L%.1f.dat",argv[1],z,global_eff,global_N_smooth,global_L/global_hubble);
      if((fid = fopen(fname,"rb"))==NULL) {
	printf("Error opening file:%s\n",fname);
	exit(1);
      }
      fread(temp,sizeof(float),global_N3_smooth,fid);  
      fclose(fid);
      for(i=0;i<global_N3_smooth;i++) {
	t21[i]=t21[i]*(1.-(double)temp[i]);  // neutral fraction...
      }
      
      sprintf(fname,"%s/deltaTb/deltaTb_z%.3lf_N%ld_L%.1f.dat",argv[1],z,global_N_smooth,global_L/global_hubble);
      if((fid = fopen(fname,"wb"))==NULL) {
	printf("Cannot open file:%s...\n",fname);
	exit(1);
      }
      for(i=0;i<global_N3_smooth;i++) temp[i]=(float)t21[i];
      fwrite(temp,sizeof(float),global_N3_smooth,fid);
      fclose(fid);
      aver=0.;
      for(i=0; i<global_N3_smooth; i++) aver+=t21[i];
      fprintf(fidtxt,"%f %.8E\n",z,aver/global_N3_smooth);
      if(z>(global_Zminsfr-global_Dzsim/10) && global_use_Lya_xrays==1) fprintf(fidtxt2,"%f %.8E\n",z,TS/global_N3_smooth);
    }
  } /* ends z cycle */
    
  fclose(fidtxt);
  exit(0);

}