/*=========================================================================== * create a timeline for current cosmological model *===========================================================================*/ void create_timeline2(double a_init, double a_fin, tlptr timeline) { int iloop; double a,t,omega,lambda,rhoc,hubble; for(iloop=0; iloop<MAXTIME2; iloop++) { a = ((double)iloop+1.0)/(double)MAXTIME2 * (a_fin-a_init) + a_init; t = calc_t(a); omega = calc_omega(a); lambda = calc_lambda(a); hubble = H0 * sqrt( simu.lambda0 * (1.-1./pow2(a)) + simu.omega0 * (1./pow3(a)-1./pow2(a)) + 1./pow2(a)); rhoc = rhoc0 * ( simu.lambda0*(1.-1./pow2(a)) + simu.omega0*(1./pow3(a)-1./pow2(a)) + 1./pow2(a)); timeline->a[iloop] = a; timeline->t[iloop] = t; timeline->omega[iloop] = omega; timeline->lambda[iloop] = lambda; timeline->hubble[iloop] = hubble; timeline->age[iloop] = t*Mpc/H0/Gyr; timeline->virial[iloop] = calc_virial(a); } }
void ahf_gridinfo(gridls *grid_list, int curgrid_no) { gridls *cur_grid; pqptr cur_pquad; cqptr cur_cquad, icur_cquad; nqptr cur_nquad, icur_nquad; nptr cur_node; int refinecounter; long x,y,z; /* The colouring tools */ int colCounter=0; int colREPLACE; int currentCol; nptr tsc_nodes[3][3][3]; /* nodes to assign to */ int firstCOLOUR; SPATIALREF *spatialRefHead; SPATIALREF *tmpSpatialRef; SPATIALREF *current; #ifdef AHFgridinfofile FILE *gridinfofile; char filename1[100]; #endif int i,tmp; int boundCount; int numWrongNodes, colRight; #ifdef AHFDEBUG char filename2[100]; FILE *fout; #endif intXYZ tmpPeriodic; int colour; SPATIALREF *previous; int iterate, startAHFgrid; double a, a3, omega, ovlim, rho_vir; double fl1dim, refine_len, refine_vol, refine_ovdens, refine_mass; double cur_shift; #ifdef VERBOSE /***************************************************************************/ fprintf(stderr,"################## ahf_gridinfo ###################\n"); fprintf(stderr,"Number of grids = %d\n", curgrid_no+1); #endif fprintf(io.logfile,"################## ahf_gridinfo ###################\n"); fprintf(io.logfile,"Number of grids = %d\n", curgrid_no+1); fflush(io.logfile); /* total number of grids to consider for AHF */ /* ahf.no_grids = curgrid_no - global.domgrid_no + 1;*/ ahf.no_grids = curgrid_no - global.domgrid_no; #ifdef VERBOSE fprintf(stderr,"Number of refinements = %d\n", ahf.no_grids); fprintf(stderr,"global.domgrid_no = %d\n", global.domgrid_no); #endif fprintf(io.logfile,"Number of refinements = %d\n", ahf.no_grids); fprintf(io.logfile,"global.domgrid_no = %d\n", global.domgrid_no); #ifdef DEBUG_GRIDS for(cur_grid=global.dom_grid+1; cur_grid<=global.dom_grid+ahf.no_grids; cur_grid++) { fprintf(stderr,"writing l1dim = %ld ... ",cur_grid->l1dim); //output_grid(cur_grid,0); write_density(cur_grid,"cur_grid-"); fprintf(stderr,"done\n"); } exit(0); #endif /* cosmology related stuff */ a = global.a; a3 = pow3(a); omega = calc_omega(a); ovlim = calc_virial(a); rho_vir = a3 * calc_rho_vir(a); /* comoving(!) density used to normalize densities */ /* get the number of the grid satisfying virial overdensity criterion */ refine_mass = simu.pmass*simu.med_weight; cur_grid = global.dom_grid; startAHFgrid = 1; for ( i=0; i<=ahf.no_grids; i++) { fl1dim = ((double)(cur_grid->l1dim)); cur_grid++; refine_len = (double)(simu.boxsize/((double)(fl1dim))); refine_vol = pow3(refine_len); refine_ovdens = (simu.Nth_ref*refine_mass/refine_vol) / rho_vir; #ifdef VERBOSE fprintf(stderr,"l1dim = %16.0f refine_ovdens = %16.4g ovlim = %16.4g\n", fl1dim,refine_ovdens,ovlim); #endif fprintf(io.logfile,"l1dim = %16.0f refine_ovdens = %16.4g ovlim = %16.4g\n", fl1dim,refine_ovdens,ovlim); if ( refine_ovdens < ovlim ) startAHFgrid = i; } /* store maximum overdensity possible with current refinement hierarchy (used by ahf_halos()!) */ global.max_ovdens = refine_ovdens; #ifdef VERBOSE fprintf(stderr,"max_ovdens = %g\n",global.max_ovdens); #endif fprintf(io.logfile,"max_ovdens = %g\n",global.max_ovdens); /* the first refinement grid to consider */ ahf.min_ref = startAHFgrid+AHF_MIN_REF_OFFSET; /* total number of grids to be used with AHF */ ahf.no_grids = ahf.no_grids - ahf.min_ref + 1; #ifdef VERBOSE fprintf(stderr,"min_ref = %d (ahf_nogrids = %d)\n", ahf.min_ref,ahf.no_grids); #endif fprintf(io.logfile,"min_ref = %d (ahf_nogrids = %d)\n", ahf.min_ref,ahf.no_grids); fflush(io.logfile); /* Initialising important variables */ #ifdef VERBOSE fprintf(stderr," initialising force.color on all grids...\n"); #endif spatialRefHead = NULL; /* Setting the HEAD of the link list of SPATIALREF's */ spatialRefTail = NULL; colREPLACE = FALSE; firstCOLOUR = 1; #ifdef WITH_OPENMP // // this loop over all grids and all refinements generates a linked-list starting // with spatialRefHead that contains the spatially connected refinements independent // of the actual grid (that information is contained in spatialRefHead->refLevel) // // to parallelize one could do this generation independently on each grid, // but would need to merge all linked-lists afterwards // //#pragma omp parallel private(cur_grid, ipquad, cur_pquad, cur_cquad, icur_cquad, cur_nquad, icur_nquad, firstCOLOUR, spatialRefHead, spatialRefTail, cur_node) shared(refinecounter) //#pragma omp for schedule(static) #endif for( refinecounter=0; refinecounter<ahf.no_grids; refinecounter++ ) { cur_grid = global.dom_grid+ahf.min_ref+refinecounter; #ifdef VERBOSE /* fprintf() to double-check that we set min_ref and ahf.no_grids correctly! */ fprintf(stderr,"%12d || l1dim = %12ld\n", refinecounter, cur_grid->l1dim); #endif /************************************************************ * Initialising all the colour tages in preperation for the * Spatial refinement identification. * Therefore, we are zeroing all the colour tags on the current refinement ************************************************************/ for(cur_pquad=cur_grid->pquad; cur_pquad != NULL; cur_pquad=cur_pquad->next) { for(cur_cquad = cur_pquad->loc; cur_cquad < cur_pquad->loc + cur_pquad->length; cur_cquad++) { for(icur_cquad = cur_cquad; icur_cquad != NULL; icur_cquad = icur_cquad->next) { for(cur_nquad = icur_cquad->loc; cur_nquad < icur_cquad->loc + icur_cquad->length; cur_nquad++) { for(icur_nquad = cur_nquad; icur_nquad != NULL; icur_nquad = icur_nquad->next) { for(cur_node = icur_nquad->loc; cur_node < icur_nquad->loc + icur_nquad->length; cur_node++) { /* zeroing the colour */ cur_node->force.colour = 0; } } } } } } /************************************************************ * Start looping through the refinement identifying the spatial refinements ************************************************************/ for(cur_pquad=cur_grid->pquad; cur_pquad != NULL; cur_pquad=cur_pquad->next) { z = cur_pquad->z; for(cur_cquad = cur_pquad->loc; cur_cquad < cur_pquad->loc + cur_pquad->length; cur_cquad++, z++) { for(icur_cquad = cur_cquad; icur_cquad != NULL; icur_cquad = icur_cquad->next) { y = icur_cquad->y; for(cur_nquad = icur_cquad->loc; cur_nquad < icur_cquad->loc + icur_cquad->length; cur_nquad++, y++) { for(icur_nquad = cur_nquad; icur_nquad != NULL; icur_nquad = icur_nquad->next) { x = icur_nquad->x; for(cur_node = icur_nquad->loc; cur_node < icur_nquad->loc + icur_nquad->length; cur_node++, x++) { /* The current colour of the node */ currentCol = cur_node->force.colour; /* We have finished updating the colours */ if ( (colREPLACE==TRUE) && (currentCol==0) ) colREPLACE=FALSE; /* Replacing the colours */ if ( colREPLACE ) { if ( col.numReplace == 2 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[1]; } else if ( col.numReplace == 3 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[0]; } else if ( col.numReplace == 4 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[0]; if (currentCol==col.holder[1]) cur_node->force.colour=col.holder[0]; } else { fprintf(stderr,"Something is wrong with colREPLACE"); } } /* New Node to check */ else { /* Getting the relevant nodes and filling colHolder */ /* search for pointers if any are null return NULL */ tsc_nodes[1][1][1] = cur_node; get_TSCnodes(cur_grid, cur_pquad, icur_cquad, icur_nquad, tsc_nodes, &z, &y, &x); /* Gathering the colour information for this node */ colourInfo(tsc_nodes); /************************************************/ /* all nodes are zero, numReplace=0 */ if ( col.numReplace == 0 ) { /* Start a new colour and keep going */ if ( (tmpSpatialRef = malloc(sizeof(SPATIALREF))) == NULL ) { fprintf(stderr,"No memory for SPATIALREF\n"); exit(1); } /* Colour the node with a unique name */ colCounter++; cur_node->force.colour = colCounter; /* Point SPATIALREF to the right information */ tmpSpatialRef->name = colCounter; tmpSpatialRef->periodic.x = 0; tmpSpatialRef->periodic.y = 0; tmpSpatialRef->periodic.z = 0; tmpSpatialRef->refLevel = refinecounter; tmpSpatialRef->cur_pquad = cur_pquad; tmpSpatialRef->cur_cquad = cur_cquad; tmpSpatialRef->icur_cquad = icur_cquad; tmpSpatialRef->cur_nquad = cur_nquad; tmpSpatialRef->icur_nquad = icur_nquad; tmpSpatialRef->cur_node = cur_node; tmpSpatialRef->x = x; tmpSpatialRef->y = y; tmpSpatialRef->z = z; /* Join the new SPATIALREF to the link list */ spatialRefTail = insertColour(tmpSpatialRef, firstCOLOUR); /* Now just for the first time we need to make head also point to this newColour */ if ( firstCOLOUR==1 ) { spatialRefHead = spatialRefTail; firstCOLOUR=0; } } /************************************************/ /* One node is non-zero, numReplace=1 */ if ( col.numReplace == 1 ) { /* set using the current colour and keep going */ cur_node->force.colour = col.holder[2]; } /************************************************/ /* Two nodes are non-zero and are not equal, numReplace=2 */ if ( col.numReplace == 2 ) { /* Replace the previous colours with the lowest colour */ /* set colour using the lowest value : in this case 1 */ cur_node->force.colour = col.holder[1]; if (!colREPLACE) { /* Find and destroy the redundant SPATIALREF */ tmpSpatialRef=deleteColour(spatialRefHead, col.holder[2]); /* STU :: This is where you start looking for the deleted colours */ cur_pquad = tmpSpatialRef->cur_pquad; cur_cquad = tmpSpatialRef->cur_cquad; icur_cquad = tmpSpatialRef->icur_cquad; cur_nquad = tmpSpatialRef->cur_nquad; icur_nquad = tmpSpatialRef->icur_nquad; cur_node = tmpSpatialRef->cur_node; x = tmpSpatialRef-> x; y = tmpSpatialRef-> y; z = tmpSpatialRef-> z; free(tmpSpatialRef); tmpSpatialRef = NULL; colREPLACE=TRUE; /* When going back to replace the colour the for loop progress the node by one - which we don't want to happen Hence, we have to do an inital sweep on the first node with the colour that is going to be deleted*/ /* Getting the relevant nodes and filling colHolder */ /* search for pointers if any are null return NULL */ currentCol = cur_node->force.colour; tsc_nodes[1][1][1] = cur_node; get_TSCnodes(cur_grid, cur_pquad, icur_cquad, icur_nquad, tsc_nodes, &z, &y, &x); /* When going back to replace the colour the for loop progress the node by one - which we don't want to happen */ if ( col.numReplace == 2 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[1]; } else if ( col.numReplace == 3 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[0]; } else if ( col.numReplace == 4 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[0]; if (currentCol==col.holder[1]) cur_node->force.colour=col.holder[0]; } else { fprintf(stderr,"Something is wrong with colREPLACE"); } /*****************************************************************/ /*****************************************************************/ } } /************************************************/ /* Two nodes are non-zero and are not equal, numReplace=2 */ /* With the last two colour elements equal */ if ( col.numReplace == 3 ) { /* Replace the previous colours with the lowest colour */ /* set colour using the lowest value : in this case 1 */ cur_node->force.colour = col.holder[0]; if (!colREPLACE) { /* Find and destroy the redendent SPATIALREF */ tmpSpatialRef=deleteColour(spatialRefHead, col.holder[2]); cur_pquad = tmpSpatialRef->cur_pquad; cur_cquad = tmpSpatialRef->cur_cquad; icur_cquad = tmpSpatialRef->icur_cquad; cur_nquad = tmpSpatialRef->cur_nquad; icur_nquad = tmpSpatialRef->icur_nquad; cur_node = tmpSpatialRef->cur_node; x = tmpSpatialRef-> x; y = tmpSpatialRef-> y; z = tmpSpatialRef-> z; free(tmpSpatialRef); tmpSpatialRef = NULL; colREPLACE=TRUE; /*****************************************************************/ /* Getting the relevant nodes and filling colHolder */ /* search for pointers if any are null return NULL */ currentCol = cur_node->force.colour; tsc_nodes[1][1][1] = cur_node; get_TSCnodes(cur_grid, cur_pquad, icur_cquad, icur_nquad, tsc_nodes, &z, &y, &x); /* When going back to replace the colour the for loop progress the node by one - which we don't want to happen */ if ( col.numReplace == 2 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[1]; } else if ( col.numReplace == 3 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[0]; } else if ( col.numReplace == 4 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[0]; if (currentCol==col.holder[1]) cur_node->force.colour=col.holder[0]; } else { fprintf(stderr,"Something is wrong with colREPLACE"); } /*****************************************************************/ /*****************************************************************/ } } /************************************************/ /* Three node are non-zero and are not equal, numReplace=3*/ if ( col.numReplace == 4 ) { /* Loop back and replace the dupulacates */ /* set colour using the lowest value : in this case 0 */ cur_node->force.colour = col.holder[0]; if (!colREPLACE) { /* Find and destroy the redendent SPATIALREF's */ tmpSpatialRef=deleteColour(spatialRefHead, col.holder[2]); free(tmpSpatialRef); tmpSpatialRef=deleteColour(spatialRefHead, col.holder[1]); cur_pquad = tmpSpatialRef->cur_pquad; cur_cquad = tmpSpatialRef->cur_cquad; icur_cquad = tmpSpatialRef->icur_cquad; cur_nquad = tmpSpatialRef->cur_nquad; icur_nquad = tmpSpatialRef->icur_nquad; cur_node = tmpSpatialRef->cur_node; x = tmpSpatialRef-> x; y = tmpSpatialRef-> y; z = tmpSpatialRef-> z; free(tmpSpatialRef); colREPLACE=TRUE; /*****************************************************************/ /*****************************************************************/ /*****************************************************************/ /*****************************************************************/ /* Getting the relevant nodes and filling colHolder */ /* search for pointers if any are null return NULL */ currentCol = cur_node->force.colour; tsc_nodes[1][1][1] = cur_node; get_TSCnodes(cur_grid, cur_pquad, icur_cquad, icur_nquad, tsc_nodes, &z, &y, &x); /* When going back to replace the colour the for loop progress the node by one - which we don't want to happen */ if ( col.numReplace == 2 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[1]; } else if ( col.numReplace == 4 ) { if (currentCol==col.holder[2]) cur_node->force.colour=col.holder[0]; if (currentCol==col.holder[1]) cur_node->force.colour=col.holder[0]; } else { fprintf(stderr,"Something is wrong with colREPLACE"); } /*****************************************************************/ /*****************************************************************/ } } /* ERROR */ if ( col.numReplace > 4 ) fprintf(stderr,"ERROR with counting col.numReplace"); } } } } } } } } /*******************************************************************************************************/ /* Checking that the grid is coloured correctly */ /*******************************************************************************************************/ #ifdef AHFDEBUG2 /* Printing the PVIEW file */ sprintf(filename2,"%shalo_pview3",global_io.params->outfile_prefix); fprintf(stderr,"%s\n",filename2); if((fout = fopen(filename2,"w")) == NULL) { fprintf(stderr,"could not open %s\n", filename2); exit(1); } #endif refinecounter=0; for( iterate=0, cur_grid=global.dom_grid+ahf.min_ref; iterate<ahf.no_grids; iterate++, cur_grid++) { /* shift of cell centre as compared to edge of box [grid units] */ cur_shift = 0.5/(double)cur_grid->l1dim; /************************************************************ * Dealing with the extents of the periodic refinements ************************************************************/ numWrongNodes = 0; boundCount = 0; for(cur_pquad=cur_grid->pquad; cur_pquad != NULL; cur_pquad=cur_pquad->next) { z = cur_pquad->z; for(cur_cquad = cur_pquad->loc; cur_cquad < cur_pquad->loc + cur_pquad->length; cur_cquad++, z++) { for(icur_cquad = cur_cquad; icur_cquad != NULL; icur_cquad = icur_cquad->next) { y = icur_cquad->y; for(cur_nquad = icur_cquad->loc; cur_nquad < icur_cquad->loc + icur_cquad->length; cur_nquad++, y++) { for(icur_nquad = cur_nquad; icur_nquad != NULL; icur_nquad = icur_nquad->next) { x = icur_nquad->x; for(cur_node = icur_nquad->loc; cur_node < icur_nquad->loc + icur_nquad->length; cur_node++, x++) { /* current node */ colour = cur_node->force.colour; tsc_nodes[1][1][1] = cur_node; get_TSCnodes(cur_grid, cur_pquad, icur_cquad, icur_nquad, tsc_nodes, &z, &y, &x); colRight = 0; colRight = checkColourInfo(tsc_nodes); numWrongNodes += colRight; #ifdef AHFDEBUG2 if ( colRight == 1 ) fprintf(fout,"p %g %g %g 0.0 1.0 0.5\n", ((float) x/(float)cur_grid->l1dim + cur_shift)*simu.boxsize, ((float) y/(float)cur_grid->l1dim + cur_shift)*simu.boxsize, ((float) z/(float)cur_grid->l1dim + cur_shift)*simu.boxsize); #endif /********************************************************************************/ /* PERIODIC BOUNDARY CONDITIONS */ /* We are on the 'bottom' boundary */ if ( (x == 0) || (y == 0) || (z == 0) ) { /* Gathering periodic information */ tmpPeriodic = testBound(tsc_nodes,x,y,z); /* This isolated refinement is also periodic */ if ( (tmpPeriodic.x == 1) || (tmpPeriodic.y == 1) || (tmpPeriodic.z == 1) ) { /* finding the corresponding isolated spatial refinement in the spatalREF linklist */ current=spatialRefHead; while ( current != NULL ) { previous = current; if ( current->name == colour ) break; current = current->next; } /* Setting the periodic elements */ if ( tmpPeriodic.x == 1 ) previous->periodic.x = 1; if ( tmpPeriodic.y == 1 ) previous->periodic.y = 1; if ( tmpPeriodic.z == 1 ) previous->periodic.z = 1; } boundCount++; } }/* current node */ } } } } } #ifdef AHFDEBUG fprintf(stderr,"Num Boundary (%d) = %d\n",refinecounter,boundCount); #endif refinecounter++; } #ifdef AHFDEBUG2 fclose(fout); #endif /**************************************************************** **************************************************************** * Creating the spatial Refinement Index * * And boundary condition information!!! * */ spatialRefIndex = NULL; if ((spatialRefIndex = calloc(colCounter+1,sizeof(SRINDEX)))==NULL) { fprintf(stderr,"Error in allocating the memory for spatialRefIndex array\n"); exit(0); } #ifdef VERBOSE fprintf(stderr,"colCounter(%d)\n",colCounter); #endif /* Creating a tempary counting array */ numIsoRef=NULL; if ((numIsoRef = calloc(ahf.no_grids,sizeof(int)))==NULL) { fprintf(stderr,"Error in allocating the memory for numIsoRef array\n"); exit(0); } for ( i=0; i<ahf.no_grids; i++) numIsoRef[i]=0; /* Setting the spatialRefIndex to -1 */ for ( i=0; i<colCounter; i++) { spatialRefIndex[i].refLevel = -1; spatialRefIndex[i].isoRefIndex = -1; } /* Loop through the spatial refinements */ /* Record their refinement level and isolated refinement index */ totnumIsoRef=0; current=spatialRefHead; while ( current!=NULL ) { spatialRefIndex[current->name].refLevel = current->refLevel; spatialRefIndex[current->name].isoRefIndex = numIsoRef[current->refLevel]; spatialRefIndex[current->name].periodic.x = current->periodic.x; spatialRefIndex[current->name].periodic.y = current->periodic.y; spatialRefIndex[current->name].periodic.z = current->periodic.z; numIsoRef[current->refLevel]++; totnumIsoRef++; current = current->next; } /**************************************************************** **************************************************************** * Looping throught the Grid structure and dumping the gridInfo file * positions - density - spatial tag */ #ifdef AHFgridinfofile fprintf(stderr,"### PRINTING GRIDINFO FILES\n"); /* Note:: do this for each refinement level starting with the first refinment */ for( iterate=0, cur_grid=global.dom_grid+ahf.min_ref; iterate<ahf.no_grids; iterate++, cur_grid++) { /* prepare filename */ sprintf(filename1,"%sz%.3f.grid.%06d",global_io.params->outfile_prefix,global.z,cur_grid->l1dim); fprintf(stderr,"%s\n",filename1); if((gridinfofile = fopen(filename1,"w")) == NULL) { fprintf(stderr,"could not open %s\n", filename1); exit(1); } /* shift of cell centre as compared to edge of box [grid units] */ cur_shift = 0.5/(double)cur_grid->l1dim; for(cur_pquad=cur_grid->pquad; cur_pquad != NULL; cur_pquad=cur_pquad->next) { z = cur_pquad->z; for(cur_cquad = cur_pquad->loc; cur_cquad < cur_pquad->loc + cur_pquad->length; cur_cquad++, z++) { for(icur_cquad = cur_cquad; icur_cquad != NULL; icur_cquad = icur_cquad->next) { y = icur_cquad->y; for(cur_nquad = icur_cquad->loc; cur_nquad < icur_cquad->loc + icur_cquad->length; cur_nquad++, y++) { for(icur_nquad = cur_nquad; icur_nquad != NULL; icur_nquad = icur_nquad->next) { x = icur_nquad->x; for(cur_node = icur_nquad->loc; cur_node < icur_nquad->loc + icur_nquad->length; cur_node++, x++) { fprintf(gridinfofile,"p %lf %lf %lf 0.0 0.0 1.0 %lf %d\n", simu.boxsize*1000*((double) x/(double)cur_grid->l1dim + cur_shift), simu.boxsize*1000*((double) y/(double)cur_grid->l1dim + cur_shift), simu.boxsize*1000*((double) z/(double)cur_grid->l1dim + cur_shift), cur_node->dens, cur_node->force.colour); } } } } } } /* Flusing and closing this grids file */ fflush(gridinfofile); fclose(gridinfofile); } #endif /* fprintf(stderr,"I'm finished here, my work is done.... time to exit :)\n"); * exit(1); */ fprintf(io.logfile,"################## ahf_gridinfo finished ##################\n\n"); fflush(io.logfile); /* free() all spatialRefHead memory */ current = spatialRefHead; while(current != NULL) { previous = current->next; free(current); current = previous; } }
/*=========================================================================== * create a timeline for current cosmological model *===========================================================================*/ void create_timeline(double a_init, double a_fin, tlptr timeline) { int iloop, iwrite_cosmo; double a,t,omega,lambda,rhoc,hubble,super_t,virial,growth,age,dummy,zred; FILE *fpcosmo; char cosmo_file[MAXSTRING], string[MAXSTRING]; #ifdef COSMOLOGY if(global_io.params->outfile_prefix != NULL) { strcpy(cosmo_file, global_io.params->outfile_prefix); strcat(cosmo_file,"cosmology"); /*--------------------------------------------- * check if there is already a .cosmology file *---------------------------------------------*/ if( (fpcosmo = fopen(cosmo_file,"r")) == NULL ) { iwrite_cosmo = TRUE; } else { iwrite_cosmo = FALSE; fclose(fpcosmo); } } else iwrite_cosmo = FALSE; #else /*--------------------------------------------- * do not write the .cosmology file whatsoever *---------------------------------------------*/ iwrite_cosmo = FALSE; #endif /* COSMOLOGY */ /*------------------------------ * create timeline from scratch *------------------------------*/ #ifndef WITH_MPI if(io.logfile != NULL) { fprintf(io.logfile,"\ncreating timeline from a=%g to a=%g (MAXTIME=%d) ... ",a_init,a_fin,MAXTIME); fflush(io.logfile); } #endif #ifdef WITH_OPENMP #pragma omp parallel private(iloop, a, t, super_t, omega, lambda, hubble, rhoc) shared(timeline, simu, a_fin, a_init) #pragma omp for schedule(static) #endif for(iloop=0; iloop<MAXTIME; iloop++) { a = ((double)iloop)/(double)(MAXTIME-1) * (a_fin-a_init) + a_init; t = calc_t(a); super_t = calc_super_t(a); omega = calc_omega(a); lambda = calc_lambda(a); hubble = calc_Hubble(a); rhoc = rhoc0 * ( simu.lambda0*(1.-1./pow2(a)) + simu.omega0*(1./pow3(a)-1./pow2(a)) + 1./pow2(a)); timeline->a[iloop] = a; timeline->t[iloop] = t; timeline->super_t[iloop] = super_t; timeline->omega[iloop] = omega; timeline->lambda[iloop] = lambda; timeline->hubble[iloop] = hubble; timeline->rhoc[iloop] = rhoc; timeline->age[iloop] = t * simu.t_unit*Mpc/Gyr; timeline->growth[iloop] = calc_growth(a); // do >>not<< call calc_virial() here when using OpenMP as calc_virial() abuses simu.omega0!!! } if(iwrite_cosmo) { fpcosmo = fopen(cosmo_file,"w"); fprintf(fpcosmo,"# z(1) a(2) t[h^-1 Gyr](3) Omega(4) lambda(5) hubble(6) RhoCrit(7) virial(8) growth(9) q(10) super_t(11)\n"); for(iloop=0; iloop<MAXTIME; iloop++) { /* calc_virial() temporarily changes simu.omega0 */ timeline->virial[iloop] = calc_virial(timeline->a[iloop]); fprintf(fpcosmo,"%10.4f %16.10f %12.6f %12.6f %12.6f %12.4f %12.6g %12.6g %12.6g %12.6g %12.6g\n", 1.0/timeline->a[iloop]-1, timeline->a[iloop], timeline->t[iloop] * simu.t_unit*Mpc/Gyr, timeline->omega[iloop], timeline->lambda[iloop], timeline->hubble[iloop], timeline->rhoc[iloop], timeline->virial[iloop], timeline->growth[iloop], calc_q(timeline->a[iloop]), timeline->super_t[iloop]); fflush(fpcosmo); } fclose(fpcosmo); } #ifndef WITH_MPI if(io.logfile != NULL) { fprintf(io.logfile,"done\n\n"); fflush(io.logfile); } #endif #ifdef DEBUG_SUPERCOMOVING { double a_prev, t_prev, super_t_prev, a_mid; double dtc, dTc; iloop = 0; a_prev = ((double)iloop)/(double)(MAXTIME-1) * (a_fin-a_init) + a_init; t_prev = calc_t(a_prev); super_t_prev = calc_super_t(a_prev); for(iloop=1; iloop<MAXTIME; iloop++) { a = ((double)iloop)/(double)(MAXTIME-1) * (a_fin-a_init) + a_init; t = calc_t(a); super_t = calc_super_t(a); dtc = t - t_prev; dTc = super_t - super_t_prev; a_mid = (a+a_prev)/2; fprintf(stderr,"z=%12.4f a=%20.10f t=%20.10f super_t=%20.10f dtc=%20.10f dTc=%20.10f dTc(dtc)=%20.10f a=%20.10f err=%g\n", 1./a-1.,a,t,super_t,dtc,dTc,dtc/(pow2(a_mid)),calc_super_a(super_t),fabs(dTc-dtc/(pow2(a_mid)))); a_prev = a; t_prev = t; super_t_prev = super_t; } } #endif #ifdef COSMOLOGY_TERM exit(0); #endif }
void do_force(FILE *fplog,t_commrec *cr, t_inputrec *inputrec, int step,t_nrnb *nrnb,gmx_wallcycle_t wcycle, gmx_localtop_t *top, gmx_groups_t *groups, matrix box,rvec x[],history_t *hist, rvec f[],rvec buf[], tensor vir_force, t_mdatoms *mdatoms, gmx_enerdata_t *enerd,t_fcdata *fcd, real lambda,t_graph *graph, t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot, real t,FILE *field,gmx_edsam_t ed, int flags) { static rvec box_size; int cg0,cg1,i,j; int start,homenr; static double mu[2*DIM]; rvec mu_tot_AB[2]; bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS,bDoForces; matrix boxs; real e,v,dvdl; t_pbc pbc; float cycles_ppdpme,cycles_pme,cycles_force; start = mdatoms->start; homenr = mdatoms->homenr; bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog)); clear_mat(vir_force); if (PARTDECOMP(cr)) { pd_cg_range(cr,&cg0,&cg1); } else { cg0 = 0; if (DOMAINDECOMP(cr)) cg1 = cr->dd->ncg_tot; else cg1 = top->cgs.nr; if (fr->n_tpi > 0) cg1--; } bStateChanged = (flags & GMX_FORCE_STATECHANGED); bNS = (flags & GMX_FORCE_NS); bFillGrid = (bNS && bStateChanged); bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr)); bDoForces = (flags & GMX_FORCE_FORCES); if (bStateChanged) { update_forcerec(fplog,fr,box); /* Calculate total (local) dipole moment in a temporary common array. * This makes it possible to sum them over nodes faster. */ calc_mu(start,homenr, x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed, mu,mu+DIM); } if (fr->ePBC != epbcNONE) { /* Compute shift vectors every step, * because of pressure coupling or box deformation! */ if (DYNAMIC_BOX(*inputrec) && bStateChanged) calc_shifts(box,fr->shift_vec); if (bCalcCGCM) { put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box, &(top->cgs),x,fr->cg_cm); inc_nrnb(nrnb,eNR_CGCM,homenr); inc_nrnb(nrnb,eNR_RESETX,cg1-cg0); } else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) { unshift_self(graph,box,x); } } else if (bCalcCGCM) { calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm); inc_nrnb(nrnb,eNR_CGCM,homenr); } if (bCalcCGCM) { if (PAR(cr)) { move_cgcm(fplog,cr,fr->cg_cm); } if (gmx_debug_at) pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr); } #ifdef GMX_MPI if (!(cr->duty & DUTY_PME)) { /* Send particle coordinates to the pme nodes. * Since this is only implemented for domain decomposition * and domain decomposition does not use the graph, * we do not need to worry about shifting. */ wallcycle_start(wcycle,ewcPP_PMESENDX); GMX_MPE_LOG(ev_send_coordinates_start); bBS = (inputrec->nwall == 2); if (bBS) { copy_mat(box,boxs); svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]); } gmx_pme_send_x(cr,bBS ? boxs : box,x,mdatoms->nChargePerturbed,lambda); GMX_MPE_LOG(ev_send_coordinates_finish); wallcycle_stop(wcycle,ewcPP_PMESENDX); } #endif /* GMX_MPI */ /* Communicate coordinates and sum dipole if necessary */ if (PAR(cr)) { wallcycle_start(wcycle,ewcMOVEX); if (DOMAINDECOMP(cr)) { dd_move_x(cr->dd,box,x,buf); } else { move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb); } /* When we don't need the total dipole we sum it in global_stat */ if (NEED_MUTOT(*inputrec)) gmx_sumd(2*DIM,mu,cr); wallcycle_stop(wcycle,ewcMOVEX); } for(i=0; i<2; i++) for(j=0;j<DIM;j++) mu_tot_AB[i][j] = mu[i*DIM + j]; if (fr->efep == efepNO) copy_rvec(mu_tot_AB[0],mu_tot); else for(j=0; j<DIM; j++) mu_tot[j] = (1.0 - lambda)*mu_tot_AB[0][j] + lambda*mu_tot_AB[1][j]; /* Reset energies */ reset_energies(&(inputrec->opts),fr,bNS,enerd,MASTER(cr)); if (bNS) { wallcycle_start(wcycle,ewcNS); if (graph && bStateChanged) /* Calculate intramolecular shift vectors to make molecules whole */ mk_mshift(fplog,graph,fr->ePBC,box,x); /* Reset long range forces if necessary */ if (fr->bTwinRange) { clear_rvecs(fr->f_twin_n,fr->f_twin); clear_rvecs(SHIFTS,fr->fshift_twin); } /* Do the actual neighbour searching and if twin range electrostatics * also do the calculation of long range forces and energies. */ dvdl = 0; ns(fplog,fr,x,f,box,groups,&(inputrec->opts),top,mdatoms, cr,nrnb,step,lambda,&dvdl,&enerd->grpp,bFillGrid,bDoForces); if (bSepDVDL) fprintf(fplog,sepdvdlformat,"LR non-bonded",0,dvdl); enerd->dvdl_lr = dvdl; enerd->term[F_DVDL] += dvdl; wallcycle_stop(wcycle,ewcNS); } if (DOMAINDECOMP(cr)) { if (!(cr->duty & DUTY_PME)) { wallcycle_start(wcycle,ewcPPDURINGPME); dd_force_flop_start(cr->dd,nrnb); } } /* Start the force cycle counter. * This counter is stopped in do_forcelow_level. * No parallel communication should occur while this counter is running, * since that will interfere with the dynamic load balancing. */ wallcycle_start(wcycle,ewcFORCE); if (bDoForces) { /* Reset PME/Ewald forces if necessary */ if (fr->bF_NoVirSum) { GMX_BARRIER(cr->mpi_comm_mygroup); if (fr->bDomDec) clear_rvecs(fr->f_novirsum_n,fr->f_novirsum); else clear_rvecs(homenr,fr->f_novirsum+start); GMX_BARRIER(cr->mpi_comm_mygroup); } /* Copy long range forces into normal buffers */ if (fr->bTwinRange) { for(i=0; i<fr->f_twin_n; i++) copy_rvec(fr->f_twin[i],f[i]); for(i=0; i<SHIFTS; i++) copy_rvec(fr->fshift_twin[i],fr->fshift[i]); } else { if (DOMAINDECOMP(cr)) clear_rvecs(cr->dd->nat_tot,f); else clear_rvecs(mdatoms->nr,f); clear_rvecs(SHIFTS,fr->fshift); } clear_rvec(fr->vir_diag_posres); GMX_BARRIER(cr->mpi_comm_mygroup); } if (inputrec->ePull == epullCONSTRAINT) clear_pull_forces(inputrec->pull); /* update QMMMrec, if necessary */ if(fr->bQMMM) update_QMMMrec(cr,fr,x,mdatoms,box,top); if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0) { /* Position restraints always require full pbc */ set_pbc(&pbc,inputrec->ePBC,box); v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms, top->idef.iparams_posres, (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres, inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda,&dvdl, fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB); if (bSepDVDL) { fprintf(fplog,sepdvdlformat, interaction_function[F_POSRES].longname,v,dvdl); } enerd->term[F_POSRES] += v; enerd->term[F_DVDL] += dvdl; inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2); } /* Compute the bonded and non-bonded forces */ do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef), cr,nrnb,wcycle,mdatoms,&(inputrec->opts), x,hist,f,enerd,fcd,box,lambda,graph,&(top->excls),mu_tot_AB, flags,&cycles_force); GMX_BARRIER(cr->mpi_comm_mygroup); if (ed) { do_flood(fplog,cr,x,f,ed,box,step); } if (DOMAINDECOMP(cr)) { dd_force_flop_stop(cr->dd,nrnb); if (wcycle) dd_cycles_add(cr->dd,cycles_force,ddCyclF); } if (bDoForces) { /* Compute forces due to electric field */ calc_f_el(MASTER(cr) ? field : NULL, start,homenr,mdatoms->chargeA,x,f,inputrec->ex,inputrec->et,t); /* When using PME/Ewald we compute the long range virial there. * otherwise we do it based on long range forces from twin range * cut-off based calculation (or not at all). */ /* Communicate the forces */ if (PAR(cr)) { wallcycle_start(wcycle,ewcMOVEF); if (DOMAINDECOMP(cr)) { dd_move_f(cr->dd,f,buf,fr->fshift); /* Position restraint do not introduce inter-cg forces */ if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl) dd_move_f(cr->dd,fr->f_novirsum,buf,NULL); } else { move_f(fplog,cr,GMX_LEFT,GMX_RIGHT,f,buf,nrnb); } wallcycle_stop(wcycle,ewcMOVEF); } } if (bDoForces) { if (vsite) { wallcycle_start(wcycle,ewcVSITESPREAD); spread_vsite_f(fplog,vsite,x,f,fr->fshift,nrnb, &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr); wallcycle_stop(wcycle,ewcVSITESPREAD); } /* Calculation of the virial must be done after vsites! */ calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f, vir_force,graph,box,nrnb,fr,inputrec->ePBC); } if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F) { /* Calculate the center of mass forces, this requires communication, * which is why pull_potential is called close to other communication. * The virial contribution is calculated directly, * which is why we call pull_potential after calc_virial. */ set_pbc(&pbc,inputrec->ePBC,box); dvdl = 0; enerd->term[F_COM_PULL] = pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc, cr,t,lambda,x,f,vir_force,&dvdl); if (bSepDVDL) fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl); enerd->term[F_DVDL] += dvdl; } if (!(cr->duty & DUTY_PME)) { cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME); dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME); } #ifdef GMX_MPI if (PAR(cr) && !(cr->duty & DUTY_PME)) { /* In case of node-splitting, the PP nodes receive the long-range * forces, virial and energy from the PME nodes here. */ wallcycle_start(wcycle,ewcPP_PMEWAITRECVF); dvdl = 0; gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl, &cycles_pme); if (bSepDVDL) fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl); enerd->term[F_COUL_RECIP] += e; enerd->term[F_DVDL] += dvdl; if (wcycle) dd_cycles_add(cr->dd,cycles_pme,ddCyclPME); wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF); } #endif if (bDoForces && fr->bF_NoVirSum) { if (vsite) { /* Spread the mesh force on virtual sites to the other particles... * This is parallellized. MPI communication is performed * if the constructing atoms aren't local. */ wallcycle_start(wcycle,ewcVSITESPREAD); spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,nrnb, &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr); wallcycle_stop(wcycle,ewcVSITESPREAD); } /* Now add the forces, this is local */ if (fr->bDomDec) { sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum); } else { sum_forces(start,start+homenr,f,fr->f_novirsum); } if (EEL_FULL(fr->eeltype)) { /* Add the mesh contribution to the virial */ m_add(vir_force,fr->vir_el_recip,vir_force); } if (debug) pr_rvecs(debug,0,"vir_force",vir_force,DIM); } /* Sum the potential energy terms from group contributions */ sum_epot(&(inputrec->opts),enerd); if (fr->print_force >= 0 && bDoForces) print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f); }