def getOrderParamq(subPos, Pos, BoxDims, lowCut=0.0, highCut=8.0): """Finds angles for 4 nearest neighbors of each water and returns for all waters the tetrahedral order parameter, q, used by Errington and Debenedetti (2001). Inputs: subPos - positions of set of atoms to measure tetrahedrality of (may be different, subset, or same as Pos) Pos - positions of ALL atoms that can make tetrahedral configurations (needed if subPos not same as Pos) BoxDims - current box dimensions to account for periodicity lowCut - lower cutoff for nearest-neighbor shell (default 0.0) highCut - higher cutoff for nearest-neighbor shell used to find 4 nearest neighbors Outputs: qVals - returns an order parameter value for each water distNeighbs - returns distances from central oxygen to 4 nearest neighbors """ #Set-up array to hold results qVals = np.zeros(len(subPos)) distNeighbs = np.zeros((len(subPos), 4)) #Find nearest neighbors for ALL atoms in subPos #But make sure using efficient algorithm... #If subPos is same as Pos, use allnearneighbors instead if np.array_equal(subPos, Pos): nearNeighbs = wl.allnearneighbors(Pos, BoxDims, lowCut, highCut).astype(bool) else: nearNeighbs = wl.nearneighbors(subPos, Pos, BoxDims, lowCut, highCut).astype(bool) #Loop over each position in subPos, finding angle made with the closest 4 neighbors, then q for (i, apos) in enumerate(subPos): #Make sure have nearest neighbors... if np.sum(nearNeighbs[i]) > 0: thisPos = wl.reimage(Pos[nearNeighbs[i]], apos, BoxDims) thisDists = np.linalg.norm(thisPos - apos, axis=1) sortInds = np.argsort(thisDists) newPos = thisPos[sortInds][:4] distNeighbs[i, :] = thisDists[sortInds][:4] #below returns symmetric, square array (zero diagonal) tempAng = wl.tetracosang(apos, newPos, BoxDims) #Only want half of array, flattened angVals = tempAng[np.triu_indices(len(tempAng), k=1)] #Now compute q for this set of angles qVals[i] = 1.0 - (3.0 / 8.0) * np.sum( (np.cos(angVals * np.pi / 180.0) + (1.0 / 3.0))**2) #Return all of the order parameter values return qVals, distNeighbs
def getCosAngs(subPos, Pos, BoxDims, lowCut=0.0, highCut=3.413): """This is called getCosAngs, but actually just returns the angles themselves (faster to convert from cos(theta) to theta in Fortran) Inputs: subPos - positions of set of atoms to measure tetrahedrality of (may be different, subset, or same as Pos) Pos - positions of ALL atoms that can make tetrahedral configurations (needed if subPos not same as Pos) BoxDims - current box dimensions to account for periodicity lowCut - lower cutoff for nearest-neighbor shell (default 0.0) highCut - higher cutoff for nearest-neighbor shell (default 3.413 - see Chaimovich, 2014, but should really change to reflect first peak in g(r) for the chosen water model) Outputs: angVals - all angle values for current configuration of positions supplied numAngs - number of angles for each central oxygen atom (i.e. number neighbors factorial) This is useful for finding which angles belong to which central oxygens This return value was added on 07/09/2017, so any code using this function before then will break, unfortunately, but the fix is easy. """ #Set-up array to hold angle results and stack as go... list increases in size! angVals = np.array([]) numAngs = np.zeros(len(subPos)) #Find nearest neighbors for ALL atoms in subPos #But make sure using efficient algorithm... #If subPos is same as Pos, use allnearneighbors instead if np.array_equal(subPos, Pos): nearNeighbs = wl.allnearneighbors(Pos, BoxDims, lowCut, highCut).astype(bool) else: nearNeighbs = wl.nearneighbors(subPos, Pos, BoxDims, lowCut, highCut).astype(bool) #Loop over each position in subPos, finding angle made with all neighbor pairs for (i, apos) in enumerate(subPos): #Make sure have nearest neighbors... if len(Pos[nearNeighbs[i]]) > 0: #below returns symmetric, square array (zero diagonal) tempAng = wl.tetracosang(apos, Pos[nearNeighbs[i]], BoxDims) #Only want half of array, flattened angVals = np.hstack( (angVals, tempAng[np.triu_indices(len(tempAng), k=1)].tolist())) numAngs[i] = tempAng.shape[0] return angVals, numAngs
def doSimDynamics(top, systemRef, integratorRef, platform, prop, temperature, scalexy=False, inBulk=False, state=None, pos=None, vels=None, nSteps=10000000): #Input a topology object, reference system, integrator, platform, platform properties, #and optionally state file, positions, or velocities #If state is specified including positions and velocities and pos and vels are not None, the #positions and velocities from the provided state will be overwritten #Does NPT, stopping periodically to run NVE to compute dynamics #Only the NPT simulation will be saved, not the NVE #Copy the reference system and integrator objects system = copy.deepcopy(systemRef) integrator = copy.deepcopy(integratorRef) #For NPT, add the barostat as a force #If not in bulk, use anisotropic barostat if not inBulk: system.addForce( mm.MonteCarloAnisotropicBarostat( (1.0, 1.0, 1.0) * u.bar, temperature, #Temperature should be SAME as for thermostat scalexy, #Set with flag for flexibility scalexy, True, #Only scale in z-direction 250 #Time-steps between MC moves )) #If in bulk, have to use isotropic barostat to avoid any weird effects with box changing dimensions else: system.addForce(mm.MonteCarloBarostat(1.0 * u.bar, temperature, 250)) #Create new simulation object for NPT simulation sim = app.Simulation(top.topology, system, integrator, platform, prop, state) #Also create copies and simulation object for the NVE we will be running systemNVE = copy.deepcopy(systemRef) integratorNVE = mm.VerletIntegrator(2.0 * u.femtoseconds) integratorNVE.setConstraintTolerance(1.0E-08) simNVE = app.Simulation(top.topology, systemNVE, integratorNVE, platform, prop) #Set the particle positions in the NPT simulation if pos is not None: sim.context.setPositions(pos) #Apply constraints before starting the simulation sim.context.applyConstraints(1.0E-08) #Check starting energy decomposition if want #decompEnergy(sim.system, sim.context.getState(getPositions=True)) #Initialize velocities if not specified if vels is not None: sim.context.setVelocities(vels) else: try: testvel = sim.context.getState(getVelocities=True).getVelocities() print( "Velocities included in state, starting with 1st particle: %s" % str(testvel[0])) #If all the velocities are zero, then set them to the temperature if not np.any(testvel.value_in_unit(u.nanometer / u.picosecond)): print( "Had velocities, but they were all zero, so setting based on temperature." ) sim.context.setVelocitiesToTemperature(temperature) except: print("Could not find velocities, setting with temperature") sim.context.setVelocitiesToTemperature(temperature) #Set up the reporter to output energies, volume, etc. sim.reporters.append( app.StateDataReporter( 'prod_out.txt', #Where to write - can be stdout or file name (default .csv, I prefer .txt) 500, #Number of steps between writes step=True, #Write step number time=True, #Write simulation time potentialEnergy=True, #Write potential energy kineticEnergy=True, #Write kinetic energy totalEnergy=True, #Write total energy temperature=True, #Write temperature volume=True, #Write volume density=False, #Write density speed=True, #Estimate of simulation speed separator= ' ' #Default is comma, but can change if want (I like spaces) )) #Set up reporter for printing coordinates (trajectory) sim.reporters.append( NetCDFReporter( 'prod.nc', #File name to write trajectory to 500, #Number of steps between writes crds=True, #Write coordinates vels=True, #Write velocities frcs=False #Write forces )) #Identify solute indices and water oxygen indices soluteInds = [] owInds = [] hw1Inds = [] hw2Inds = [] for res in top.residues: if res.name not in ['OTM', 'CTM', 'STM', 'NTM', 'SOL']: for atom in res.atoms: soluteInds.append(atom.idx) elif res.name == 'SOL': for atom in res.atoms: if atom.name == 'OW': owInds.append(atom.idx) elif atom.name == 'HW1': hw1Inds.append(atom.idx) elif atom.name == 'HW2': hw2Inds.append(atom.idx) print("Solute indices:") print(soluteInds) #print("Water oxygen indices:") #print(owInds) #print("Water hydrogen (1st) indices:") #print(hw1Inds) #print("Water hydrogen (2nd) indices:") #print(hw2Inds) #Define cutoffs for solute solvation shells solShell1Cut = 0.55 #nanometers from all solute atoms (including hydrogens) solShell2Cut = 0.85 #Create array to store the dynamic information of interest every 0.2 ps (100 steps) for 50 ps calcSteps = 100 calcTotSteps = 25000 numWats = np.zeros( (int(calcTotSteps / calcSteps) + 1, 2) ) #Number waters that started in shell that are in shell at later time dipCorrs = np.zeros((int(calcTotSteps / calcSteps) + 1, 2)) #Dipole correlation in both solute shells #Start running dynamics print("\nRunning NPT simulation with interspersed NVE to find dynamics...") sim.context.setTime(0.0) stepChunk = 5000 #Run NVE for 50 ps to find dynamics every 10 ps countSteps = 0 while countSteps < nSteps: countSteps += stepChunk sim.step(stepChunk) #Record the simulation state so can kick off the NVE simulation thisState = sim.context.getState(getPositions=True, getVelocities=True) #Get solute and water oxygen coordinates after wrapping around the solute coords = thisState.getPositions(asNumpy=True) boxDims = np.diagonal(thisState.getPeriodicBoxVectors(asNumpy=True)) wrapCOM = np.average(coords[soluteInds], axis=0) coords = wl.reimage(coords, wrapCOM, boxDims) - wrapCOM solCoords = coords[soluteInds] owCoords = coords[owInds] hw1Coords = coords[hw1Inds] hw2Coords = coords[hw2Inds] #Figure out which waters are in the solute solvation shells shell1BoolMat = wl.nearneighbors(solCoords, owCoords, boxDims, 0.0, solShell1Cut) shell1Bool = np.array(np.sum(shell1BoolMat, axis=0), dtype=bool) shell2BoolMat = wl.nearneighbors(solCoords, owCoords, boxDims, solShell1Cut, solShell2Cut) shell2Bool = np.array(np.sum(shell2BoolMat, axis=0), dtype=bool) #Count number of waters in each shell (will need for averaging) thisCount1 = int(np.sum(shell1Bool)) thisCount2 = int(np.sum(shell2Bool)) #print("Found %i waters in shell1"%thisCount1) #print("Found %i waters in shell2"%thisCount2) #Loop over waters in shells and compute dipole vectors as references refDipoles1 = np.zeros((thisCount1, 3)) refDipoles2 = np.zeros((thisCount2, 3)) for k, pos in enumerate(owCoords[shell1Bool]): thisOHvecs = wl.reimage( [hw1Coords[shell1Bool][k], hw2Coords[shell1Bool][k]], pos, boxDims) - pos thisDip = -0.5 * (thisOHvecs[0] + thisOHvecs[1]) refDipoles1[k] = thisDip / np.linalg.norm(thisDip) for k, pos in enumerate(owCoords[shell2Bool]): thisOHvecs = wl.reimage( [hw1Coords[shell2Bool][k], hw2Coords[shell2Bool][k]], pos, boxDims) - pos thisDip = -0.5 * (thisOHvecs[0] + thisOHvecs[1]) refDipoles2[k] = thisDip / np.linalg.norm(thisDip) #Set up the NVE simulation simNVE.context.setState(thisState) simNVE.context.setTime(0.0) #Loop over taking steps to computed dynamics countStepsNVE = 0 while countStepsNVE <= calcTotSteps: calcState = simNVE.context.getState(getPositions=True) #Get solute and water oxygen coordinates after wrapping around the solute coords = calcState.getPositions(asNumpy=True) wrapCOM = np.average(coords[soluteInds], axis=0) coords = wl.reimage(coords, wrapCOM, boxDims) - wrapCOM solCoords = coords[soluteInds] owCoords = coords[owInds] hw1Coords = coords[hw1Inds] hw2Coords = coords[hw2Inds] #Count waters that started in each shell that are now in the shell at this time #No absorbing boundaries thisbool1Mat = wl.nearneighbors(solCoords, owCoords, boxDims, 0.0, solShell1Cut) thisbool1 = np.array(np.sum(thisbool1Mat, axis=0), dtype=bool) thisbool2Mat = wl.nearneighbors(solCoords, owCoords, boxDims, solShell1Cut, solShell2Cut) thisbool2 = np.array(np.sum(thisbool2Mat, axis=0), dtype=bool) numWats[int(countStepsNVE / calcSteps), 0] += int(np.sum(thisbool1 * shell1Bool)) numWats[int(countStepsNVE / calcSteps), 1] += int(np.sum(thisbool2 * shell2Bool)) #Loop over waters in shells and compute dipole vectors for this configuration #Adding to sum that we will normalize to find average at each time point for k, pos in enumerate(owCoords[shell1Bool]): thisOHvecs = wl.reimage( [hw1Coords[shell1Bool][k], hw2Coords[shell1Bool][k]], pos, boxDims) - pos thisDip = -0.5 * (thisOHvecs[0] + thisOHvecs[1]) thisDip /= np.linalg.norm(thisDip) dipCorrs[int(countStepsNVE / calcSteps), 0] += (np.dot(thisDip, refDipoles1[k]) / float(thisCount1)) for k, pos in enumerate(owCoords[shell2Bool]): thisOHvecs = wl.reimage( [hw1Coords[shell2Bool][k], hw2Coords[shell2Bool][k]], pos, boxDims) - pos thisDip = -0.5 * (thisOHvecs[0] + thisOHvecs[1]) thisDip /= np.linalg.norm(thisDip) dipCorrs[int(countStepsNVE / calcSteps), 1] += (np.dot(thisDip, refDipoles2[k]) / float(thisCount2)) simNVE.step(calcSteps) countStepsNVE += calcSteps #Finish normalizing dipole correlations (really cosine of angle between dipole vector at different times) numWats /= float(int(nSteps / stepChunk)) dipCorrs /= float(int(nSteps / stepChunk)) print("Normalizing factor for finding averages: %f" % float(int(nSteps / stepChunk))) #And save the final state of the NPT simulation in case we want to extend it sim.saveState('nptDynamicsState.xml') #And return the dipole correlations and times at which they were computed timeVals = 0.002 * np.arange(0.0, calcTotSteps + 0.0001, calcSteps) return numWats, dipCorrs, timeVals
def main(args): print time.ctime(time.time()) #Get topology file we're working with topFile = args[0] #And figure out if we're dealing with solute at surface or in bulk if (args[1] == 'True'): inBulk = True else: inBulk = False if (args[2] == 'True'): doReweight = True else: doReweight = False #Read in topology file now to get information on solute atoms top = pmd.load_file(topFile) soluteInds = [] for res in top.residues: if res.name not in ['OTM', 'CTM', 'STM', 'NTM', 'SOL']: for atom in res.atoms: soluteInds.append(atom.idx) #Now define how we compute three-body angles with bins and cut-off #Shell cut-off shellCut = 3.32 #1st minimum distance for TIP4P-Ew water at 298.15 K and 1 bar #Number of angle bins nAngBins = 100 #500 #Define bin centers (should be nBins equally spaced between 0 and 180) angBinCents = 0.5 * (np.arange(0.0, 180.001, 180.0/nAngBins)[:-1] + np.arange(0.0, 180.001, 180.0/nAngBins)[1:]) #And distance bins for local oxygen-oxygen RDF calculation #(really distance histograms from central oxygens - can normalize however we want, really) distBinWidth = 0.05 nDistBins = int(shellCut / distBinWidth) distBins = np.arange(0.0, nDistBins*distBinWidth+0.00001, distBinWidth) distBinCents = 0.5 * (distBins[:-1] + distBins[1:]) #Define the size of the probes used for assessing density fluctuations near solute probeRadius = 3.3 # radius in Angstroms; the DIAMETER of a methane so assumes other atoms methane-sized #And bins for numbers of waters in probes probeBins = np.arange(0.0, 21.00001, 1.0) nProbeBins = len(probeBins) - 1 #Will use np.histogram, which includes left edge in bin (so if want up to 20, go to 21) #And will record number waters in each solvation shell (histograms of) shellBins = np.arange(0.0, 251.00001, 1.0) #probably way too many bins, but don't want to run out nShellBins = len(shellBins) - 1 #Should we do 2D histogram for angles and distances? #Or also do 2D histogram for number of waters in probe and three-body angle? #Interesting if make probe radius same size as three-body angle cutoff? #Finally, also define the bins for computing RDFs of waters near all solute atoms #Will use to define 1st and 2nd solvation shells rdfBinWidth = 0.2 rdfMax = 12.00 rdfBins = np.arange(0.0, rdfMax+0.000001, rdfBinWidth) rdfBinCents = 0.5 * (rdfBins[:-1] + rdfBins[1:]) nRDFBins = len(rdfBinCents) rdfBinVols = (4.0*np.pi/3.0)*(rdfBins[1:]**3 - rdfBins[:-1]**3) bulkDens = 0.0332 #Roughly right for TIP4P-EW at 298.15 K and 1 bar in inverse Angstroms cubed #Need to create a variety of arrays to hold the data we're interested in #Will record distributions in both 1st and 2nd solvation shells shellCountsCoupled = np.zeros((nShellBins, 2)) #Histograms for numbers of waters in hydration shells of solutes probeHistsCoupled = np.zeros((nProbeBins, 2)) #Histograms for numbers waters in probes in 1st and 2nd hydration shells angHistsCoupled = np.zeros((nAngBins, 2)) #Histograms of three-body angles for water oxygens within solvation shells distHistsCoupled = np.zeros((nDistBins, 2)) #Histograms of distances to water oxygens from central oxygens solRDFsCoupled = np.zeros((nRDFBins, len(soluteInds))) #RDFs between each solute atom and water oxygens shellCountsDecoupled = np.zeros((nShellBins, 2)) #Same as above, but decoupled state, not coupled state probeHistsDecoupled = np.zeros((nProbeBins, 2)) angHistsDecoupled = np.zeros((nAngBins, 2)) distHistsDecoupled = np.zeros((nDistBins, 2)) solRDFsDecoupled = np.zeros((nRDFBins, len(soluteInds))) #First need configuration weights to use in computing average quantities #But only do if we're using a simuation with an expanded ensemble if doReweight: if inBulk: weightsCoupled, weightsDecoupled = getConfigWeightsBulk(kB=1.0, T=1.0) #Using 1 for kB and T because alchemical_output.txt should already have potential energies in kBT simDirs = ['.'] else: weightsCoupled, weightsDecoupled = getConfigWeightsSurf() simDirs = ['Quad_0.25X_0.25Y', 'Quad_0.25X_0.75Y', 'Quad_0.75X_0.25Y', 'Quad_0.75X_0.75Y'] else: weightsCoupled = np.array([]) weightsDecoupled = np.array([]) simDirs = ['.'] #To correctly match weights up to configurations, need to count frames from all trajectories countFrames = 0 #Next, want to loop over all trajectories and compute RDFs from solute atoms to water oxygens #Will use this to define solvation shells for finding other properties #Actually, having looked at RDFs, just use 5.5 for first shell and 8.5 for second shell... #AND use all atoms, including hydrogens, which have LJ interactions in GAFF2, to define shells #Actually now only using heavy atoms... but when look at RDFs, examine all atoms for adir in simDirs: if doReweight: #Before loading trajectory, figure out how many frames to exclude due to weight equilibration alcDat = np.loadtxt(adir+'/alchemical_output.txt') startTime = alcDat[0, 1] startFrame = int(startTime) - 1 else: startFrame = 0 top = pmd.load_file(topFile) top.rb_torsions = pmd.TrackedList([]) #This is just for SAM systems so that it doesn't break pytraj top = pt.load_parmed(top, traj=False) traj = pt.iterload(adir+'/prod.nc', top, frame_slice=(startFrame, -1)) if not doReweight: weightsCoupled = np.hstack((weightsCoupled, np.ones(len(traj)))) weightsDecoupled = np.hstack((weightsDecoupled, np.ones(len(traj)))) print("\nTopology and trajectory loaded from directory %s" % adir) owInds = top.select('@OW') soluteInds = top.select('!(:OTM,CTM,STM,NTM,SOL)') print("\n\tFound %i water oxygens" % len(owInds)) print("\tFound %i solute atoms" % len(soluteInds)) for i, frame in enumerate(traj): if i%1000 == 0: print "On frame %i" % i boxDims = np.array(frame.box.values[:3]) currCoords = np.array(frame.xyz) #Wrap based on soluate atom center of geometry and get coordinates of interest wrapCOM = np.average(currCoords[soluteInds], axis=0) currCoords = wl.reimage(currCoords, wrapCOM, boxDims) - wrapCOM owCoords = currCoords[owInds] solCoords = currCoords[soluteInds] #Loop over solute atoms and find pair-distance histograms with water oxygens for j, acoord in enumerate(solCoords): solRDFsCoupled[:,j] += (weightsCoupled[countFrames+i] * wl.pairdistancehistogram(np.array([acoord]), owCoords, rdfBinWidth, nRDFBins, boxDims)) solRDFsDecoupled[:,j] += (weightsDecoupled[countFrames+i] * wl.pairdistancehistogram(np.array([acoord]), owCoords, rdfBinWidth, nRDFBins, boxDims)) #Note that pairdistancehistogram is right-edge inclusive, NOT left-edge inclusive #In practice, not a big difference countFrames += len(traj) #Finish by normalizing RDFs properly for j in range(len(soluteInds)): solRDFsCoupled[:,j] /= rdfBinVols #bulkDens*rdfBinVols solRDFsDecoupled[:,j] /= rdfBinVols #bulkDens*rdfBinVols if not doReweight: solRDFsCoupled /= float(countFrames) solRDFsDecoupled /= float(countFrames) #And save to file np.savetxt('solute-OW_RDFs_coupled.txt', np.hstack((np.array([rdfBinCents]).T, solRDFsCoupled)), header='RDF bins (A) solute atom-OW RDF for solute atom indices %s'%(str(soluteInds))) np.savetxt('solute-OW_RDFs_decoupled.txt', np.hstack((np.array([rdfBinCents]).T, solRDFsDecoupled)), header='RDF bins (A) solute atom-OW RDF for solute atom indices %s'%(str(soluteInds))) print("\tFound RDFs for water oxygens from solute indices.") solShell1Cut = 5.5 #Angstroms from all solute atoms (including hydrogens) solShell2Cut = 8.5 #And now that we know how many frames, we can assign real weights if not reweighting if not doReweight: weightsCoupled /= float(countFrames) weightsDecoupled /= float(countFrames) #Reset countFrames so get weights right countFrames = 0 #Repeat looping over trajectories to calculate water properties in solute solvation shell for adir in simDirs: if doReweight: #Before loading trajectory, figure out how many frames to exclude due to weight equilibration alcDat = np.loadtxt(adir+'/alchemical_output.txt') startTime = alcDat[0, 1] startFrame = int(startTime) - 1 else: startFrame = 0 top = pmd.load_file(topFile) top.rb_torsions = pmd.TrackedList([]) #This is just for SAM systems so that it doesn't break pytraj top = pt.load_parmed(top, traj=False) traj = pt.iterload(adir+'/prod.nc', top, frame_slice=(startFrame, -1)) print("\nTopology and trajectory loaded from directory %s" % adir) owInds = top.select('@OW') soluteInds = top.select('!(:OTM,CTM,STM,NTM,SOL)&!(@H=)') surfInds = top.select('(:OTM,CTM,STM,NTM)&!(@H=)') #For probe insertions, also include solute and surface heavy atoms print("\n\tFound %i water oxygens" % len(owInds)) print("\tFound %i solute heavy atoms" % len(soluteInds)) print("\tFound %i non-hydrogen surface atoms" % len(surfInds)) if len(surfInds) == 0: surfInds.dtype=int for i, frame in enumerate(traj): #if i%10 == 0: # print "On frame %i" % i boxDims = np.array(frame.box.values[:3]) currCoords = np.array(frame.xyz) #Wrap based on soluate atom center of geometry and get coordinates of interest wrapCOM = np.average(currCoords[soluteInds], axis=0) currCoords = wl.reimage(currCoords, wrapCOM, boxDims) - wrapCOM owCoords = currCoords[owInds] solCoords = currCoords[soluteInds] surfCoords = currCoords[surfInds] #Now get solvent shells around solute shell1BoolMat = wl.nearneighbors(solCoords, owCoords, boxDims, 0.0, solShell1Cut) shell1Bool = np.array(np.sum(shell1BoolMat, axis=0), dtype=bool) shell2BoolMat = wl.nearneighbors(solCoords, owCoords, boxDims, solShell1Cut, solShell2Cut) shell2Bool = np.array(np.sum(shell2BoolMat, axis=0), dtype=bool) #And add weight to histogram for numbers of waters in shells thisCount1 = int(np.sum(shell1Bool)) shellCountsCoupled[thisCount1, 0] += weightsCoupled[countFrames+i] shellCountsDecoupled[thisCount1, 0] += weightsDecoupled[countFrames+i] thisCount2 = int(np.sum(shell2Bool)) shellCountsCoupled[thisCount2, 1] += weightsCoupled[countFrames+i] shellCountsDecoupled[thisCount2, 1] += weightsDecoupled[countFrames+i] #And compute water properties of solvent shells, first 3-body angles thisAngs1, thisNumAngs1 = wp.getCosAngs(owCoords[shell1Bool], owCoords, boxDims, highCut=shellCut) thisAngHist1, thisAngBins1 = np.histogram(thisAngs1, bins=nAngBins, range=[0.0, 180.0], density=False) angHistsCoupled[:,0] += weightsCoupled[countFrames+i] * thisAngHist1 angHistsDecoupled[:,0] += weightsDecoupled[countFrames+i] * thisAngHist1 thisAngs2, thisNumAngs2 = wp.getCosAngs(owCoords[shell2Bool], owCoords, boxDims, highCut=shellCut) thisAngHist2, thisAngBins2 = np.histogram(thisAngs2, bins=nAngBins, range=[0.0, 180.0], density=False) angHistsCoupled[:,1] += weightsCoupled[countFrames+i] * thisAngHist2 angHistsDecoupled[:,1] += weightsDecoupled[countFrames+i] * thisAngHist2 #And ow-ow pair distance histograms in both shells as well thisDistHist1 = wl.pairdistancehistogram(owCoords[shell1Bool], owCoords, distBinWidth, nDistBins, boxDims) distHistsCoupled[:,0] += weightsCoupled[countFrames+i] * thisDistHist1 distHistsDecoupled[:,0] += weightsDecoupled[countFrames+i] * thisDistHist1 thisDistHist2 = wl.pairdistancehistogram(owCoords[shell2Bool], owCoords, distBinWidth, nDistBins, boxDims) distHistsCoupled[:,1] += weightsCoupled[countFrames+i] * thisDistHist2 distHistsDecoupled[:,1] += weightsDecoupled[countFrames+i] * thisDistHist2 #Next compute distributions of numbers of waters in probes centered within each shell #To do this, create random grid of points in SQUARE that encompasses both shells #Then only keep points within each shell based on distance #Square will be based on shell cutoffs and min and max coordinates in each dimension of solute minSolX = np.min(solCoords[:,0]) - solShell2Cut maxSolX = np.max(solCoords[:,0]) + solShell2Cut minSolY = np.min(solCoords[:,1]) - solShell2Cut maxSolY = np.max(solCoords[:,1]) + solShell2Cut minSolZ = np.min(solCoords[:,2]) - solShell2Cut maxSolZ = np.max(solCoords[:,2]) + solShell2Cut thisGridX = minSolX + np.random.random(500)*(maxSolX - minSolX) thisGridY = minSolY + np.random.random(500)*(maxSolY - minSolY) thisGridZ = minSolZ + np.random.random(500)*(maxSolZ - minSolZ) thisGrid = np.vstack((thisGridX, thisGridY, thisGridZ)).T gridBoolMat1 = wl.nearneighbors(solCoords, thisGrid, boxDims, 0.0, solShell1Cut) gridBool1 = np.array(np.sum(gridBoolMat1, axis=0), dtype=bool) thisNum1 = wl.probegrid(np.vstack((owCoords, surfCoords, solCoords)), thisGrid[gridBool1], probeRadius, boxDims) thisProbeHist1, thisProbeBins1 = np.histogram(thisNum1, bins=probeBins, density=False) probeHistsCoupled[:,0] += weightsCoupled[countFrames+i] * thisProbeHist1 probeHistsDecoupled[:,0] += weightsDecoupled[countFrames+i] * thisProbeHist1 gridBoolMat2 = wl.nearneighbors(solCoords, thisGrid, boxDims, solShell1Cut, solShell2Cut) gridBool2 = np.array(np.sum(gridBoolMat2, axis=0), dtype=bool) thisNum2 = wl.probegrid(np.vstack((owCoords, surfCoords, solCoords)), thisGrid[gridBool2], probeRadius, boxDims) thisProbeHist2, thisProbeBins2 = np.histogram(thisNum2, bins=probeBins, density=False) probeHistsCoupled[:,1] += weightsCoupled[countFrames+i] * thisProbeHist2 probeHistsDecoupled[:,1] += weightsDecoupled[countFrames+i] * thisProbeHist2 countFrames += len(traj) #Should have everything we need, so save to text files np.savetxt('solute_shell_hists.txt', np.hstack((np.array([shellBins[:-1]]).T, shellCountsCoupled, shellCountsDecoupled)), header='Histograms of numbers of waters in first and second solute solvation shells with solvent in coupled (columns 2, 3) and decoupled (columns 4, 5) states') np.savetxt('solute_probe_hists.txt', np.hstack((np.array([probeBins[:-1]]).T, probeHistsCoupled, probeHistsDecoupled)), header='Number waters in probe histograms in first and second solute solvation shells with solvent in coupled (columns 2, 3) and decoupled (columns 4, 5) states') np.savetxt('solute_ang_hists.txt', np.hstack((np.array([angBinCents]).T, angHistsCoupled, angHistsDecoupled)), header='3-body angle histograms in first and second solute solvation shells with solvent in coupled (columns 2, 3) and decoupled (columns 4, 5) states') np.savetxt('solute_pair_hists.txt', np.hstack((np.array([distBinCents]).T, distHistsCoupled, distHistsDecoupled)), header='O-O pair-distance histograms in first and second solute solvation shells with solvent in coupled (columns 2, 3) and decoupled (columns 4, 5) states') print time.ctime(time.time())
def computeSphericalFourierCoeffs(subPos, Pos, BoxDims, lowCut=0.0, highCut=3.413, minDegree=0, maxDegree=12): """Computes the vectors of Fourier coefficients for each degree of a spherical harmonic expansion as described by Keys, Iacovella, and Glotzer, 2011. subPos is treated as the central atoms, while Pos should include the atoms that may potentially be neighbors. Inputs: subPos - positions of atoms to treat as the central atoms Pos - positions of all other atoms, which will be considered for neighbor-searching; can be same as subPos BoxDims - box dimensions so that minimum images may be used lowCut - (default 0.0) the lower cutoff for the radial shell highCut - (default 3.413) the upper cutoff for the radial shell minDegree - (default 0) the minimum spherical harmonic degree (l) maxDegree - (default 12) the maximum spherical harmonic degree (l) Outputs: coeffVecs - a len(subPos) x (1 + maxDegree - minDegree) x (2*maxDegree + 1) matrix For each central atom in subPos, a matrix of the complex-valued vectors (as rows) is provided. This still allows magnitudes to be easily evaluated, since real and imaginary parts of zero will contribute nothing to the magnitude numNeighbs - number of neighbors for each water molecule (necessary to compute global order parameters or coefficients by multiplying by this and the dividing by the total of the waters to "average" over) """ #Set up the output matrix now since know size coeffVecs = np.zeros( (len(subPos), 1 + maxDegree - minDegree, 2 * maxDegree + 1), dtype=complex) #And array to return number of neighbors for each water numNeighbs = np.zeros(len(subPos), dtype='float16') #Would be nice to combine neighbor searching with 3-body angle computation or H-bonding code #But then harder to efficiently implement different cutoffs... #So that might be too ambitious #Find neighbors within cutoff for ALL atoms in subPos #But make sure using efficient algorithm... #If subPos is same as Pos, use allnearneighbors instead if np.array_equal(subPos, Pos): nearNeighbs = wl.allnearneighbors(Pos, BoxDims, lowCut, highCut).astype(bool) else: nearNeighbs = wl.nearneighbors(subPos, Pos, BoxDims, lowCut, highCut).astype(bool) #Loop over each position in subPos and find neighbor positions in spherical coordinates for (i, apos) in enumerate(subPos): #Make sure have nearest neighbors... if len(Pos[nearNeighbs[i]]) > 0: tempPos = wl.reimage(Pos[nearNeighbs[i]], apos, BoxDims) - apos numNeighbs[i] = len(tempPos) #Compute radial distances... unfortunate that have to do this again, but maybe improve later rdists = np.linalg.norm(tempPos, axis=1) #And get polar and azimuthal angles polarang = np.arccos(tempPos[:, 2] / rdists) azimang = np.arctan2( tempPos[:, 1], tempPos[:, 0] ) #Using special arctan2 function to get quadrant right #Now compute Fourier coefficient vectors (i.e. have complex-valued component of coefficient vector #associated with each m value, where m = -l, -l+1, ... , l) #Loop over the desired number of coefficients to compute for l in range(minDegree, maxDegree + 1): thisvec = np.zeros(2 * l + 1, dtype=complex) #Also note that have one vector for each neighbor, so must loop over neighbors for j in range(len(tempPos)): thisvec += sph_harm(np.arange(-l, l + 1), l, azimang[j], polarang[j]) thisvec /= len(tempPos) #And compute the magnitude of this vector of complex numbers coeffVecs[i, l - minDegree, :(2 * l + 1)] = thisvec return coeffVecs, numNeighbs
def main(args): print('\nReading in files and obtaining LJ information...') #First read in topology and trajectory topFile = args[0] trajFile = args[1] top = pmd.load_file(topFile) trajtop = copy.deepcopy(top) trajtop.rb_torsions = pmd.TrackedList( []) #Necessary for SAM systems so doesn't break pytraj trajtop = pt.load_parmed(trajtop, traj=False) traj = pt.iterload(trajFile, top=trajtop) #Next use parmed to get LJ parameters for all atoms in the solute, as well as water oxygens and surface atoms #While go, also collect dictionaries of atomic indices associated with each atom type #Will have to check separately when looking at overlaps with solute soluteLJ = [] soluteInds = [] dictOtherLJ = {} dictOtherInds = {} for res in top.residues: if res.name not in ['OTM', 'CTM', 'STM', 'NTM', 'SOL']: for atom in res.atoms: soluteInds.append(atom.idx) soluteLJ.append([atom.sigma, atom.epsilon]) else: for atom in res.atoms: if not atom.type in dictOtherInds: dictOtherInds[atom.type] = [atom.idx] else: dictOtherInds[atom.type].append(atom.idx) if not atom.type in dictOtherLJ: dictOtherLJ[atom.type] = np.array( [atom.sigma, atom.epsilon]) soluteLJ = np.array(soluteLJ) print(soluteLJ) #Use Lorentz-Berthelot combining rules to get LJ parameters between each solute atom and a water oxygen dictMixLJ = {} for i in range(soluteLJ.shape[0]): for akey in dictOtherLJ: dictMixLJ['%i_%s' % (i, akey)] = np.array([ 0.5 * (soluteLJ[i, 0] + dictOtherLJ[akey][0]), np.sqrt(soluteLJ[i, 1] * dictOtherLJ[akey][1]) ]) for key, val in dictMixLJ.iteritems(): print("%s, %s" % (key, str(val.tolist()))) print( '\nDetermining hard-sphere radii for all combinations of solute and other system atoms...' ) #Next compute hard-sphere radii by integrating according to Barker and Hendersen, Weeks, etc. #In order to this right, technically using WCA potential, not LJ hsRadii = {} rvals = np.arange(0.0, 50.005, 0.005) betaLJ = 1.0 / ((1.9872036E-03) * 298.15) for i in range(soluteLJ.shape[0]): for akey in dictOtherLJ: [thisSig, thisEps] = dictMixLJ['%i_%s' % (i, akey)] if thisEps == 0.0: hsRadii['%i_%s' % (i, akey)] = 0.0 continue thisRmin = thisSig * (2.0**(1.0 / 6.0)) thisSigdRmin6 = (thisSig / thisRmin)**6 thisPotRmin = 4.0 * thisEps * ((thisSigdRmin6**2) - thisSigdRmin6) thisSigdR6 = (thisSig / rvals)**6 thisPotVals = 4.0 * thisEps * ( (thisSigdR6**2) - thisSigdR6) - thisPotRmin thisPotVals[np.where(rvals >= thisRmin)] = 0.0 thisPotVals *= betaLJ thisIntegrand = 1.0 - np.exp(-thisPotVals) thisIntegrand[0] = 1.0 hsRadii['%i_%s' % (i, akey)] = simps(thisIntegrand, rvals) #Need to multiply by two because will use atom centers to check overlap? Is Rhs distance between centers? for key, val in hsRadii.iteritems(): print("%s, %f" % (key, val)) #Keep track of hard-sphere radii with water oxygens specially solOWhsRadii = np.zeros(soluteLJ.shape[0]) for i in range(soluteLJ.shape[0]): solOWhsRadii[i] = hsRadii['%i_OW_tip4pew' % i] #Only using TIP4P/EW here print('\nStarting loop over trajectory...') #Now loop over trajectory and check if solute is overlapping with any waters OR surface atoms #Will create a distribution of overlapping atoms numOverlap = np.arange(101.0) countOverlap = np.zeros(len(numOverlap)) #And also track average solute volume that we're trying to insert solVol = 0.0 countFrames = 0 print('') for frame in traj: if countFrames % 100 == 0: print("On frame %i" % countFrames) countFrames += 1 boxDims = np.array(frame.box.values[:3]) currCoords = np.array(frame.xyz) #Get solute coordinates and make sure solute is whole #Then get solute coordinates relative to first solute atom #Don't need any other wrapping because wl.nearneighbors will do its own wrapping solCoords = currCoords[soluteInds] solRefCoords = wl.reimage(solCoords, solCoords[0], boxDims) - solCoords[0] #While we have the coordinates nice, compute solute volume #Do this based on hard-sphere radii to water oxygens, which is the most likely case anyway solVol += np.sum(wl.spherevolumes(solRefCoords, solOWhsRadii, 0.1)) #Will shift the solute first atom (and others too) to a number of random locations #Since applying restraint (if have surface) Z is drawn from the distribution we want #X and Y can be drawn more easily from uniform distributions, so do this to randomize solute position numRandXY = 1000 randX = np.random.random(numRandXY) * boxDims[0] randY = np.random.random(numRandXY) * boxDims[1] #And will keep track of number of overlapping atoms for each solute random position thisTotOverlap = np.zeros(numRandXY, dtype=int) #Loop over all non-solute atoms in the system by atom type for akey, val in dictOtherInds.iteritems(): #For this specific atom type, need to keep track of WHICH neighbors overlapping #And need to do for EACH solute atom #Don't want to double-count if two solute atoms both overlap the same other atom overlapBool = np.zeros((numRandXY, len(val)), dtype=int) #Now loop over each solute atom for i, coord in enumerate(solRefCoords): thisRadius = hsRadii['%i_%s' % (i, akey)] if thisRadius == 0.0: continue #Define coordinates of the current solute atom we're working with #So setting first atom to random XY position, then shifting by distance to first atom hsCoords = np.zeros((numRandXY, 3)) hsCoords[:, 0] = coord[0] + randX hsCoords[:, 1] = coord[1] + randY hsCoords[:, 2] = coord[2] + solCoords[0, 2] #Identify boolean for overlapping atoms and add to overall boolean for overlap #Note that want OR operation, so adding boolean arrays overlapBool += wl.nearneighbors(hsCoords, currCoords[val], boxDims, 0.0, thisRadius) #For this non-solute atom type, add number of atoms overlapping with ANY solute atom thisTotOverlap += np.sum(np.array(overlapBool, dtype=bool), axis=1) thisBins = np.arange(np.max(thisTotOverlap) + 1) countOverlap[thisBins] += np.bincount(thisTotOverlap) print(countOverlap.tolist()) print('Hard-sphere solute insertion probability: %f' % (-np.log(countOverlap[0] / np.sum(countOverlap)))) #Save the distribution to file np.savetxt( 'HS-solute_overlap_hist.txt', np.vstack((numOverlap, countOverlap)).T, header= 'Number of non-solute atoms overlapping Histogram count') solVol /= float(len(traj)) print( 'Average solute hard-sphere volume (based on water oxygen LJ params): %f' % (solVol))