int gmx_sans(int argc, char *argv[]) { const char *desc[] = { "[THISMODULE] computes SANS spectra using Debye formula.", "It currently uses topology file (since it need to assigne element for each atom).", "[PAR]", "Parameters:[PAR]" "[TT]-pr[tt] Computes normalized g(r) function averaged over trajectory[PAR]", "[TT]-prframe[tt] Computes normalized g(r) function for each frame[PAR]", "[TT]-sq[tt] Computes SANS intensity curve averaged over trajectory[PAR]", "[TT]-sqframe[tt] Computes SANS intensity curve for each frame[PAR]", "[TT]-startq[tt] Starting q value in nm[PAR]", "[TT]-endq[tt] Ending q value in nm[PAR]", "[TT]-qstep[tt] Stepping in q space[PAR]", "Note: When using Debye direct method computational cost increases as", "1/2 * N * (N - 1) where N is atom number in group of interest.", "[PAR]", "WARNING: If sq or pr specified this tool can produce large number of files! Up to two times larger than number of frames!" }; static gmx_bool bPBC = TRUE; static gmx_bool bNORM = FALSE; static real binwidth = 0.2, grid = 0.05; /* bins shouldnt be smaller then smallest bond (~0.1nm) length */ static real start_q = 0.0, end_q = 2.0, q_step = 0.01; static real mcover = -1; static unsigned int seed = 0; static int nthreads = -1; static const char *emode[] = { NULL, "direct", "mc", NULL }; static const char *emethod[] = { NULL, "debye", "fft", NULL }; gmx_neutron_atomic_structurefactors_t *gnsf; gmx_sans_t *gsans; #define NPA asize(pa) t_pargs pa[] = { { "-bin", FALSE, etREAL, {&binwidth}, "[HIDDEN]Binwidth (nm)" }, { "-mode", FALSE, etENUM, {emode}, "Mode for sans spectra calculation" }, { "-mcover", FALSE, etREAL, {&mcover}, "Monte-Carlo coverage should be -1(default) or (0,1]"}, { "-method", FALSE, etENUM, {emethod}, "[HIDDEN]Method for sans spectra calculation" }, { "-pbc", FALSE, etBOOL, {&bPBC}, "Use periodic boundary conditions for computing distances" }, { "-grid", FALSE, etREAL, {&grid}, "[HIDDEN]Grid spacing (in nm) for FFTs" }, {"-startq", FALSE, etREAL, {&start_q}, "Starting q (1/nm) "}, {"-endq", FALSE, etREAL, {&end_q}, "Ending q (1/nm)"}, { "-qstep", FALSE, etREAL, {&q_step}, "Stepping in q (1/nm)"}, { "-seed", FALSE, etINT, {&seed}, "Random seed for Monte-Carlo"}, #ifdef GMX_OPENMP { "-nt", FALSE, etINT, {&nthreads}, "Number of threads to start"}, #endif }; FILE *fp; const char *fnTPX, *fnNDX, *fnTRX, *fnDAT = NULL; t_trxstatus *status; t_topology *top = NULL; t_atom *atom = NULL; gmx_rmpbc_t gpbc = NULL; gmx_bool bTPX; gmx_bool bFFT = FALSE, bDEBYE = FALSE; gmx_bool bMC = FALSE; int ePBC = -1; matrix box; char title[STRLEN]; rvec *x; int natoms; real t; char **grpname = NULL; atom_id *index = NULL; int isize; int i, j; char *hdr = NULL; char *suffix = NULL; t_filenm *fnmdup = NULL; gmx_radial_distribution_histogram_t *prframecurrent = NULL, *pr = NULL; gmx_static_structurefactor_t *sqframecurrent = NULL, *sq = NULL; output_env_t oenv; #define NFILE asize(fnm) t_filenm fnm[] = { { efTPX, "-s", NULL, ffREAD }, { efTRX, "-f", NULL, ffREAD }, { efNDX, NULL, NULL, ffOPTRD }, { efDAT, "-d", "nsfactor", ffOPTRD }, { efXVG, "-pr", "pr", ffWRITE }, { efXVG, "-sq", "sq", ffWRITE }, { efXVG, "-prframe", "prframe", ffOPTWR }, { efXVG, "-sqframe", "sqframe", ffOPTWR } }; nthreads = gmx_omp_get_max_threads(); if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv)) { return 0; } /* check that binwidth not smaller than smallers distance */ check_binwidth(binwidth); check_mcover(mcover); /* setting number of omp threads globaly */ gmx_omp_set_num_threads(nthreads); /* Now try to parse opts for modes */ switch (emethod[0][0]) { case 'd': bDEBYE = TRUE; switch (emode[0][0]) { case 'd': bMC = FALSE; break; case 'm': bMC = TRUE; break; default: break; } break; case 'f': bFFT = TRUE; break; default: break; } if (bDEBYE) { if (bMC) { fprintf(stderr, "Using Monte Carlo Debye method to calculate spectrum\n"); } else { fprintf(stderr, "Using direct Debye method to calculate spectrum\n"); } } else if (bFFT) { gmx_fatal(FARGS, "FFT method not implemented!"); } else { gmx_fatal(FARGS, "Unknown combination for mode and method!"); } /* Try to read files */ fnDAT = ftp2fn(efDAT, NFILE, fnm); fnTPX = ftp2fn(efTPX, NFILE, fnm); fnTRX = ftp2fn(efTRX, NFILE, fnm); gnsf = gmx_neutronstructurefactors_init(fnDAT); fprintf(stderr, "Read %d atom names from %s with neutron scattering parameters\n\n", gnsf->nratoms, fnDAT); snew(top, 1); snew(grpname, 1); snew(index, 1); bTPX = read_tps_conf(fnTPX, title, top, &ePBC, &x, NULL, box, TRUE); printf("\nPlease select group for SANS spectra calculation:\n"); get_index(&(top->atoms), ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, grpname); gsans = gmx_sans_init(top, gnsf); /* Prepare reference frame */ if (bPBC) { gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr); gmx_rmpbc(gpbc, top->atoms.nr, box, x); } natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box); if (natoms != top->atoms.nr) { fprintf(stderr, "\nWARNING: number of atoms in tpx (%d) and trajectory (%d) do not match\n", natoms, top->atoms.nr); } do { if (bPBC) { gmx_rmpbc(gpbc, top->atoms.nr, box, x); } /* allocate memory for pr */ if (pr == NULL) { /* in case its first frame to read */ snew(pr, 1); } /* realy calc p(r) */ prframecurrent = calc_radial_distribution_histogram(gsans, x, box, index, isize, binwidth, bMC, bNORM, mcover, seed); /* copy prframecurrent -> pr and summ up pr->gr[i] */ /* allocate and/or resize memory for pr->gr[i] and pr->r[i] */ if (pr->gr == NULL) { /* check if we use pr->gr first time */ snew(pr->gr, prframecurrent->grn); snew(pr->r, prframecurrent->grn); } else { /* resize pr->gr and pr->r if needed to preven overruns */ if (prframecurrent->grn > pr->grn) { srenew(pr->gr, prframecurrent->grn); srenew(pr->r, prframecurrent->grn); } } pr->grn = prframecurrent->grn; pr->binwidth = prframecurrent->binwidth; /* summ up gr and fill r */ for (i = 0; i < prframecurrent->grn; i++) { pr->gr[i] += prframecurrent->gr[i]; pr->r[i] = prframecurrent->r[i]; } /* normalize histo */ normalize_probability(prframecurrent->grn, prframecurrent->gr); /* convert p(r) to sq */ sqframecurrent = convert_histogram_to_intensity_curve(prframecurrent, start_q, end_q, q_step); /* print frame data if needed */ if (opt2fn_null("-prframe", NFILE, fnm)) { snew(hdr, 25); snew(suffix, GMX_PATH_MAX); /* prepare header */ sprintf(hdr, "g(r), t = %f", t); /* prepare output filename */ fnmdup = dup_tfn(NFILE, fnm); sprintf(suffix, "-t%.2f", t); add_suffix_to_output_names(fnmdup, NFILE, suffix); fp = xvgropen(opt2fn_null("-prframe", NFILE, fnmdup), hdr, "Distance (nm)", "Probability", oenv); for (i = 0; i < prframecurrent->grn; i++) { fprintf(fp, "%10.6f%10.6f\n", prframecurrent->r[i], prframecurrent->gr[i]); } done_filenms(NFILE, fnmdup); fclose(fp); sfree(hdr); sfree(suffix); sfree(fnmdup); } if (opt2fn_null("-sqframe", NFILE, fnm)) { snew(hdr, 25); snew(suffix, GMX_PATH_MAX); /* prepare header */ sprintf(hdr, "I(q), t = %f", t); /* prepare output filename */ fnmdup = dup_tfn(NFILE, fnm); sprintf(suffix, "-t%.2f", t); add_suffix_to_output_names(fnmdup, NFILE, suffix); fp = xvgropen(opt2fn_null("-sqframe", NFILE, fnmdup), hdr, "q (nm^-1)", "s(q)/s(0)", oenv); for (i = 0; i < sqframecurrent->qn; i++) { fprintf(fp, "%10.6f%10.6f\n", sqframecurrent->q[i], sqframecurrent->s[i]); } done_filenms(NFILE, fnmdup); fclose(fp); sfree(hdr); sfree(suffix); sfree(fnmdup); } /* free pr structure */ sfree(prframecurrent->gr); sfree(prframecurrent->r); sfree(prframecurrent); /* free sq structure */ sfree(sqframecurrent->q); sfree(sqframecurrent->s); sfree(sqframecurrent); } while (read_next_x(oenv, status, &t, x, box)); close_trj(status); /* normalize histo */ normalize_probability(pr->grn, pr->gr); sq = convert_histogram_to_intensity_curve(pr, start_q, end_q, q_step); /* prepare pr.xvg */ fp = xvgropen(opt2fn_null("-pr", NFILE, fnm), "G(r)", "Distance (nm)", "Probability", oenv); for (i = 0; i < pr->grn; i++) { fprintf(fp, "%10.6f%10.6f\n", pr->r[i], pr->gr[i]); } xvgrclose(fp); /* prepare sq.xvg */ fp = xvgropen(opt2fn_null("-sq", NFILE, fnm), "I(q)", "q (nm^-1)", "s(q)/s(0)", oenv); for (i = 0; i < sq->qn; i++) { fprintf(fp, "%10.6f%10.6f\n", sq->q[i], sq->s[i]); } xvgrclose(fp); /* * Clean up memory */ sfree(pr->gr); sfree(pr->r); sfree(pr); sfree(sq->q); sfree(sq->s); sfree(sq); please_cite(stdout, "Garmay2012"); return 0; }
int gmx_genconf(int argc, char *argv[]) { const char *desc[] = { "[THISMODULE] multiplies a given coordinate file by simply stacking them", "on top of each other, like a small child playing with wooden blocks.", "The program makes a grid of [IT]user-defined[it]", "proportions ([TT]-nbox[tt]), ", "and interspaces the grid point with an extra space [TT]-dist[tt].[PAR]", "When option [TT]-rot[tt] is used the program does not check for overlap", "between molecules on grid points. It is recommended to make the box in", "the input file at least as big as the coordinates + ", "van der Waals radius.[PAR]", "If the optional trajectory file is given, conformations are not", "generated, but read from this file and translated appropriately to", "build the grid." }; const char *bugs[] = { "The program should allow for random displacement of lattice points." }; int vol; t_atoms *atoms; /* list with all atoms */ rvec *x, *xx, *v; /* coordinates? */ real t; vec4 *xrot, *vrot; int ePBC; matrix box, boxx; /* box length matrix */ rvec shift; int natoms; /* number of atoms in one molecule */ int nres; /* number of molecules? */ int i, j, k, l, m, ndx, nrdx, nx, ny, nz; t_trxstatus *status; gmx_bool bTRX; gmx_output_env_t *oenv; t_filenm fnm[] = { { efSTX, "-f", "conf", ffREAD }, { efSTO, "-o", "out", ffWRITE }, { efTRX, "-trj", NULL, ffOPTRD } }; #define NFILE asize(fnm) static rvec nrbox = {1, 1, 1}; static int seed = 0; /* seed for random number generator */ static gmx_bool bRandom = FALSE; /* False: no random rotations */ static gmx_bool bRenum = TRUE; /* renumber residues */ static rvec dist = {0, 0, 0}; /* space added between molecules ? */ static rvec max_rot = {180, 180, 180}; /* maximum rotation */ t_pargs pa[] = { { "-nbox", FALSE, etRVEC, {nrbox}, "Number of boxes" }, { "-dist", FALSE, etRVEC, {dist}, "Distance between boxes" }, { "-seed", FALSE, etINT, {&seed}, "Random generator seed (0 means generate)" }, { "-rot", FALSE, etBOOL, {&bRandom}, "Randomly rotate conformations" }, { "-maxrot", FALSE, etRVEC, {max_rot}, "Maximum random rotation" }, { "-renumber", FALSE, etBOOL, {&bRenum}, "Renumber residues" } }; if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, asize(bugs), bugs, &oenv)) { return 0; } if (seed == 0) { seed = static_cast<int>(gmx::makeRandomSeed()); } gmx::DefaultRandomEngine rng(seed); bTRX = ftp2bSet(efTRX, NFILE, fnm); nx = (int)(nrbox[XX]+0.5); ny = (int)(nrbox[YY]+0.5); nz = (int)(nrbox[ZZ]+0.5); if ((nx <= 0) || (ny <= 0) || (nz <= 0)) { gmx_fatal(FARGS, "Number of boxes (-nbox) should be larger than zero"); } vol = nx*ny*nz; /* calculate volume in grid points (= nr. molecules) */ t_topology *top; snew(top, 1); atoms = &top->atoms; read_tps_conf(opt2fn("-f", NFILE, fnm), top, &ePBC, &x, &v, box, FALSE); natoms = atoms->nr; nres = atoms->nres; /* nr of residues in one element? */ /* make space for all the atoms */ add_t_atoms(atoms, natoms*(vol-1), nres*(vol-1)); srenew(x, natoms*vol); /* get space for coordinates of all atoms */ srenew(v, natoms*vol); /* velocities. not really needed? */ snew(xrot, natoms); /* get space for rotation matrix? */ snew(vrot, natoms); if (bTRX) { if (!read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &xx, boxx)) { gmx_fatal(FARGS, "No atoms in trajectory %s", ftp2fn(efTRX, NFILE, fnm)); } } else { snew(xx, natoms); for (i = 0; i < natoms; i++) { copy_rvec(x[i], xx[i]); } } for (k = 0; (k < nz); k++) /* loop over all gridpositions */ { shift[ZZ] = k*(dist[ZZ]+box[ZZ][ZZ]); for (j = 0; (j < ny); j++) { shift[YY] = j*(dist[YY]+box[YY][YY])+k*box[ZZ][YY]; for (i = 0; (i < nx); i++) { shift[XX] = i*(dist[XX]+box[XX][XX])+j*box[YY][XX]+k*box[ZZ][XX]; ndx = (i*ny*nz+j*nz+k)*natoms; nrdx = (i*ny*nz+j*nz+k)*nres; /* Random rotation on input coords */ if (bRandom) { rand_rot(natoms, xx, v, xrot, vrot, &rng, max_rot); } for (l = 0; (l < natoms); l++) { for (m = 0; (m < DIM); m++) { if (bRandom) { x[ndx+l][m] = xrot[l][m]; v[ndx+l][m] = vrot[l][m]; } else { x[ndx+l][m] = xx[l][m]; v[ndx+l][m] = v[l][m]; } } if (ePBC == epbcSCREW && i % 2 == 1) { /* Rotate around x axis */ for (m = YY; m <= ZZ; m++) { x[ndx+l][m] = box[YY][m] + box[ZZ][m] - x[ndx+l][m]; v[ndx+l][m] = -v[ndx+l][m]; } } for (m = 0; (m < DIM); m++) { x[ndx+l][m] += shift[m]; } atoms->atom[ndx+l].resind = nrdx + atoms->atom[l].resind; atoms->atomname[ndx+l] = atoms->atomname[l]; } for (l = 0; (l < nres); l++) { atoms->resinfo[nrdx+l] = atoms->resinfo[l]; if (bRenum) { atoms->resinfo[nrdx+l].nr += nrdx; } } if (bTRX) { if (!read_next_x(oenv, status, &t, xx, boxx) && ((i+1)*(j+1)*(k+1) < vol)) { gmx_fatal(FARGS, "Not enough frames in trajectory"); } } } } } if (bTRX) { close_trj(status); } /* make box bigger */ for (m = 0; (m < DIM); m++) { box[m][m] += dist[m]; } svmul(nx, box[XX], box[XX]); svmul(ny, box[YY], box[YY]); svmul(nz, box[ZZ], box[ZZ]); if (ePBC == epbcSCREW && nx % 2 == 0) { /* With an even number of boxes in x we can forgot about the screw */ ePBC = epbcXYZ; } /*depending on how you look at it, this is either a nasty hack or the way it should work*/ if (bRenum) { for (i = 0; i < atoms->nres; i++) { atoms->resinfo[i].nr = i+1; } } write_sto_conf(opt2fn("-o", NFILE, fnm), *top->name, atoms, x, v, ePBC, box); sfree(x); sfree(v); sfree(xrot); sfree(vrot); sfree(xx); done_top(top); sfree(top); done_filenms(NFILE, fnm); output_env_done(oenv); return 0; }