def initialize_from_hdf5_file(file_name, structure, read_trajectory=True, initial_cut=1, final_cut=None, memmap=False): import h5py print("Reading data from hdf5 file: " + file_name) trajectory = None velocity = None vc = None reduced_q_vector = None #Check file exists if not os.path.isfile(file_name): print(file_name + ' file does not exist!') exit() hdf5_file = h5py.File(file_name, "r") if "trajectory" in hdf5_file and read_trajectory is True: trajectory = hdf5_file['trajectory'][:] if final_cut is not None: trajectory = trajectory[initial_cut-1:final_cut] else: trajectory = trajectory[initial_cut-1:] if "velocity" in hdf5_file: velocity = hdf5_file['velocity'][:] if final_cut is not None: velocity = velocity[initial_cut-1:final_cut] else: velocity = velocity[initial_cut-1:] if "vc" in hdf5_file: vc = hdf5_file['vc'][:] if final_cut is not None: vc = vc[initial_cut-1:final_cut] else: vc = vc[initial_cut-1:] if "reduced_q_vector" in hdf5_file: reduced_q_vector = hdf5_file['reduced_q_vector'][:] print("Load trajectory projected onto {0}".format(reduced_q_vector)) time = hdf5_file['time'][:] supercell = hdf5_file['super_cell'][:] hdf5_file.close() if vc is None: return dyn.Dynamics(structure=structure, trajectory=trajectory, velocity=velocity, time=time, supercell=np.dot(np.diagflat(supercell), structure.get_cell()), memmap=memmap) else: return vc, reduced_q_vector, dyn.Dynamics(structure=structure, time=time, supercell=np.dot(np.diagflat(supercell), structure.get_cell()), memmap=memmap)
def generate_test_trajectory(structure, supercell=(1, 1, 1), minimum_frequency=0.1, # THz total_time=2, # picoseconds time_step=0.002, # picoseconds temperature=400, # Kelvin silent=False, memmap=False, phase_0=0.0): import random from dynaphopy.power_spectrum import _progress_bar print('Generating ideal harmonic data for testing') kb_boltzmann = 0.831446 # u * A^2 / ( ps^2 * K ) number_of_unit_cells_phonopy = np.prod(np.diag(structure.get_supercell_phonon())) number_of_unit_cells = np.prod(supercell) # atoms_relation = float(number_of_unit_cells)/ number_of_unit_cells_phonopy #Recover dump trajectory from file (test only) import pickle if False: dump_file = open( "trajectory.save", "r" ) trajectory = pickle.load(dump_file) return trajectory number_of_atoms = structure.get_number_of_cell_atoms() number_of_primitive_atoms = structure.get_number_of_primitive_atoms() number_of_dimensions = structure.get_number_of_dimensions() positions = structure.get_positions(supercell=supercell) masses = structure.get_masses(supercell=supercell) number_of_atoms = number_of_atoms*number_of_unit_cells number_of_primitive_cells = number_of_atoms/number_of_primitive_atoms atom_type = structure.get_atom_type_index(supercell=supercell) #Generate additional wave vectors sample # structure.set_supercell_phonon_renormalized(np.diag(supercell)) q_vector_list = pho_interface.get_commensurate_points(structure, np.diag(supercell)) q_vector_list_cart = [ np.dot(q_vector, 2*np.pi*np.linalg.inv(structure.get_primitive_cell()).T) for q_vector in q_vector_list] atoms_relation = float(len(q_vector_list)*number_of_primitive_atoms)/number_of_atoms #Generate frequencies and eigenvectors for the testing wave vector samples print('Wave vectors included in test (commensurate points)') eigenvectors_r = [] frequencies_r = [] for i in range(len(q_vector_list)): print(q_vector_list[i]) eigenvectors, frequencies = pho_interface.obtain_eigenvectors_and_frequencies(structure, q_vector_list[i]) eigenvectors_r.append(eigenvectors) frequencies_r.append(frequencies) number_of_frequencies = len(frequencies_r[0]) #Generating trajectory if not silent: _progress_bar(0, 'generating') #Generating trajectory trajectory = [] for time in np.arange(total_time, step=time_step): coordinates = np.array(positions[:, :], dtype=complex) for i_freq in range(number_of_frequencies): for i_long, q_vector in enumerate(q_vector_list_cart): if abs(frequencies_r[i_long][i_freq]) > minimum_frequency: # Prevent error due to small frequencies amplitude = np.sqrt(number_of_dimensions * kb_boltzmann * temperature / number_of_primitive_cells * atoms_relation)/(frequencies_r[i_long][i_freq] * 2 * np.pi) # + random.uniform(-1,1)*0.05 normal_mode = amplitude * np.exp(np.complex(0, -1) * frequencies_r[i_long][i_freq] * 2.0 * np.pi * time) phase = np.exp(np.complex(0, 1) * np.dot(q_vector, positions.T) + phase_0) coordinates += (1.0 / np.sqrt(masses)[None].T * eigenvectors_r[i_long][i_freq, atom_type] * phase[None].T * normal_mode).real trajectory.append(coordinates) if not silent: _progress_bar(float(time + time_step) / total_time, 'generating', ) trajectory = np.array(trajectory) time = np.array([i * time_step for i in range(trajectory.shape[0])], dtype=float) energy = np.array([number_of_atoms * number_of_dimensions * kb_boltzmann * temperature for i in range(trajectory.shape[0])], dtype=float) #Save a trajectory object to file for later recovery (test only) if False: dump_file = open("trajectory.save", "w") pickle.dump(dyn.Dynamics(structure=structure, trajectory=np.array(trajectory, dtype=complex), energy=np.array(energy), time=time, supercell=np.dot(np.diagflat(supercell), structure.get_cell())), dump_file) dump_file.close() # structure.set_supercell_phonon_renormalized(None) return dyn.Dynamics(structure=structure, trajectory=np.array(trajectory,dtype=complex), energy=np.array(energy), time=time, supercell=np.dot(np.diagflat(supercell), structure.get_cell()), memmap=memmap)
def read_from_file_test(): print('Reading structure from test file') #Condicions del test number_of_dimensions = 2 f_coordinates = open('Data Files/test.out', 'r') f_velocity = open('Data Files/test2.out', 'r') f_trajectory = open('Data Files/test3.out', 'r') #Coordinates reading positions = [] while True: row = f_coordinates.readline().split() if not row: break for i in range(len(row)): row[i] = float(row[i]) positions.append(row) atom_type = np.array(positions,dtype=int)[:, 2] positions = np.array(positions)[:,:number_of_dimensions] print('Coordinates reading complete') structure = atomtest.Structure(positions=positions, atomic_numbers=atom_type, cell=[[2,0],[0,1]], masses=[1] * positions.shape[0]) #all 1 number_of_atoms = structure.get_number_of_atoms() structure.set_number_of_primitive_atoms(2) print('number of atoms in primitive cell') print(structure.get_number_of_primitive_atoms()) print('number of total atoms in structure (super cell)') print(number_of_atoms) #Velocity reading section velocity = [] while True: row = f_velocity.readline().replace('I','j').replace('*','').replace('^','E').split() if not row: break for i in range(len(row)): row[i] = complex('('+row[i]+')') velocity.append(row) # velocity = velocity[:4000][:] #Limitate the number of points (just for testing) time = np.array([velocity[i][0] for i in range(len(velocity))]).real velocity = np.array([[[velocity[i][j*number_of_dimensions+k+1] for k in range(number_of_dimensions)] for j in range(number_of_atoms)] for i in range (len(velocity))]) print('Velocity reading complete') #Trajectory reading trajectory = [] while True: row = f_trajectory.readline().replace('I','j').replace('*','').replace('^','E').split() if not row: break for i in range(len(row)): row[i] = complex('('+row[i]+')') trajectory.append(row) trajectory = np.array([[[trajectory[i][j*number_of_dimensions+k+1] for k in range(number_of_dimensions)] for j in range(number_of_atoms)] for i in range (len(trajectory))]) print('Trajectory reading complete') return dyn.Dynamics(trajectory=trajectory, #velocity=velocity, time=time, structure=structure)
def generate_lammps_trajectory( structure, input_file, total_time=0.1, # picoseconds time_step=0.002, # picoseconds relaxation_time=0, silent=False, supercell=(1, 1, 1), memmap=False, # not fully implemented yet! velocity_only=False, lammps_log=True, sampling_interval=1): # in timesteps cmdargs_lammps = ['-echo', 'none', '-screen', 'none'] if not lammps_log: cmdargs_lammps += ['-log', 'none'] lmp = lammps(cmdargs=cmdargs_lammps) lmp.file(input_file) lmp.command('timestep {}'.format(time_step)) lmp.command('replicate {} {} {}'.format(*supercell)) lmp.command('run 0') # natoms = lmp.extract_global("natoms",0) # mass = lmp.extract_atom("mass",2) # temp = lmp.extract_compute("thermo_temp",0,0) # print("Temperature from compute =",temp) # print("Natoms, mass, x[0][0] coord =", natoms, mass[1], x[0][0]) # print ('thermo', lmp.get_thermo('1')) xlo = lmp.extract_global("boxxlo", 1) xhi = lmp.extract_global("boxxhi", 1) ylo = lmp.extract_global("boxylo", 1) yhi = lmp.extract_global("boxyhi", 1) zlo = lmp.extract_global("boxzlo", 1) zhi = lmp.extract_global("boxzhi", 1) xy = lmp.extract_global("xy", 1) yz = lmp.extract_global("yz", 1) xz = lmp.extract_global("xz", 1) simulation_cell = np.array([[xhi - xlo, xy, xz], [0, yhi - ylo, yz], [0, 0, zhi - zlo]]).T positions = [] velocity = [] energy = [] na = lmp.get_natoms() xc = lmp.gather_atoms("x", 1, 3) reference = np.array([xc[i] for i in range(na * 3)]).reshape((na, 3)) template = get_correct_arrangement(reference, structure) indexing = np.argsort(template) lmp.command('variable energy equal pe'.format( int(relaxation_time / time_step))) lmp.command('run {}'.format(int(relaxation_time / time_step))) if not silent: _progress_bar(0, 'lammps') n_loops = int(total_time / time_step / sampling_interval) for i in range(n_loops): if not silent: _progress_bar( float((i + 1) * time_step * sampling_interval) / total_time, 'lammps', ) lmp.command('run {}'.format(sampling_interval)) xc = lmp.gather_atoms("x", 1, 3) vc = lmp.gather_atoms("v", 1, 3) energy.append(lmp.extract_variable('energy', 'all', 0)) velocity.append( np.array([vc[i] for i in range(na * 3)]).reshape( (na, 3))[indexing, :]) if not velocity_only: positions.append( np.array([xc[i] for i in range(na * 3)]).reshape( (na, 3))[indexing, :]) positions = np.array(positions, dtype=complex) velocity = np.array(velocity, dtype=complex) energy = np.array(energy) if velocity_only: positions = None lmp.close() time = np.array( [i * time_step * sampling_interval for i in range(n_loops)], dtype=float) return dyn.Dynamics(structure=structure, trajectory=positions, velocity=velocity, time=time, energy=energy, supercell=simulation_cell, memmap=memmap)
def generate_lammps_trajectory( structure, input_file, total_time=0.1, # picoseconds time_step=0.002, # picoseconds relaxation_time=0, silent=False, supercell=(1, 1, 1), memmap=False, # not fully implemented yet! velocity_only=False, lammps_log=True, temperature=None, thermostat_mass=0.5, sampling_interval=1): # in timesteps cmdargs_lammps = ['-echo', 'none', '-screen', 'none'] if not lammps_log: cmdargs_lammps += ['-log', 'none'] lmp = lammps(cmdargs=cmdargs_lammps) lmp.file(input_file) lmp.command('timestep {}'.format(time_step)) lmp.command('replicate {} {} {}'.format(*supercell)) # Force reset temperature (overwrites lammps script) # This forces NVT simulation if temperature is not None: print('Temperature reset to {} (NVT)'.format(temperature)) lmp.command('fix int all nvt temp {0} {0} {1}'.format( temperature, thermostat_mass)) # Check if correct number of atoms if lmp.extract_global("natoms", 0) < 2: print('Number of atoms in MD should be higher than 1!') exit() # Check if initial velocities all zero if not np.array(lmp.gather_atoms("v", 1, 3)).any(): print('Set lammps initial velocities') t = temperature if temperature is not None else 100 lmp.command( 'velocity all create {} 3627941 dist gaussian mom yes'. format(t)) lmp.command('velocity all scale {}'.format(t)) lmp.command('run 0') # natoms = lmp.extract_global("natoms",0) # mass = lmp.extract_atom("mass",2) # temp = lmp.extract_compute("thermo_temp",0,0) # print("Temperature from compute =",temp) # print("Natoms, mass, x[0][0] coord =", natoms, mass[1], x[0][0]) # print ('thermo', lmp.get_thermo('1')) try: xlo = lmp.extract_global("boxxlo") xhi = lmp.extract_global("boxxhi") ylo = lmp.extract_global("boxylo") yhi = lmp.extract_global("boxyhi") zlo = lmp.extract_global("boxzlo") zhi = lmp.extract_global("boxzhi") xy = lmp.extract_global("xy") yz = lmp.extract_global("yz") xz = lmp.extract_global("xz") except TypeError: xlo = lmp.extract_global("boxxlo", 1) xhi = lmp.extract_global("boxxhi", 1) ylo = lmp.extract_global("boxylo", 1) yhi = lmp.extract_global("boxyhi", 1) zlo = lmp.extract_global("boxzlo", 1) zhi = lmp.extract_global("boxzhi", 1) xy = lmp.extract_global("xy", 1) yz = lmp.extract_global("yz", 1) xz = lmp.extract_global("xz", 1) except UnboundLocalError: boxlo, boxhi, xy, yz, xz, periodicity, box_change = lmp.extract_box() xlo, ylo, zlo = boxlo xhi, yhi, zhi = boxhi simulation_cell = np.array([[xhi - xlo, xy, xz], [0, yhi - ylo, yz], [0, 0, zhi - zlo]]).T positions = [] velocity = [] energy = [] na = lmp.get_natoms() id = lmp.extract_atom("id", 0) id = np.array([id[i] - 1 for i in range(na)], dtype=int) xp = lmp.extract_atom("x", 3) reference = np.array([[xp[i][0], xp[i][1], xp[i][2]] for i in range(na)], dtype=float)[id, :] # xc = lmp.gather_atoms("x", 1, 3) # reference = np.array([xc[i] for i in range(na*3)]).reshape((na,3)) template = get_correct_arrangement(reference, structure) indexing = np.argsort(template) lmp.command('variable energy equal pe'.format( int(relaxation_time / time_step))) lmp.command('run {}'.format(int(relaxation_time / time_step))) if not silent: _progress_bar(0, 'lammps') n_loops = int(total_time / time_step / sampling_interval) for i in range(n_loops): if not silent: _progress_bar( float((i + 1) * time_step * sampling_interval) / total_time, 'lammps', ) lmp.command('run {}'.format(sampling_interval)) # xc = lmp.gather_atoms("x", 1, 3) # vc = lmp.gather_atoms("v", 1, 3) energy.append(lmp.extract_variable('energy', 'all', 0)) id = lmp.extract_atom("id", 0) id = np.array([id[i] - 1 for i in range(na)], dtype=int) vc = lmp.extract_atom("v", 3) velocity.append( np.array([[vc[i][0], vc[i][1], vc[i][2]] for i in range(na)], dtype=float)[id, :][indexing, :]) # velocity.append(np.array([vc[i] for i in range(na * 3)]).reshape((na, 3))[indexing, :]) if not velocity_only: xc = lmp.extract_atom("x", 3) positions.append( np.array([[xc[i][0], xc[i][1], xc[i][2]] for i in range(na)], dtype=float)[id, :][indexing, :]) # positions.append(np.array([xc[i] for i in range(na * 3)]).reshape((na, 3))[indexing, :]) positions = np.array(positions, dtype=complex) velocity = np.array(velocity, dtype=complex) energy = np.array(energy) if velocity_only: positions = None lmp.close() time = np.array( [i * time_step * sampling_interval for i in range(n_loops)], dtype=float) return dyn.Dynamics(structure=structure, trajectory=positions, velocity=velocity, time=time, energy=energy, supercell=simulation_cell, memmap=memmap)
def read_vasp_trajectory( file_name, structure=None, time_step=None, limit_number_steps=10000000, # Maximum number of steps read (for security) last_steps=None, initial_cut=1, end_cut=None, memmap=False, template=None): # warning warnings.warn( 'This parser will be deprecated, you can use XDATCAR instead', DeprecationWarning) # Check file exists if not os.path.isfile(file_name): print('Trajectory file does not exist!') exit() # Check time step if time_step is not None: print( 'Warning! Time step flag has no effect reading from VASP OUTCAR file (time step will be read from file)' ) if memmap: print( 'Warning! Memory mapping is not implemented in VASP OUTCAR parser') # Starting reading print("Reading VASP trajectory") print("This could take long, please wait..") # Dimensionality of VASP calculation number_of_dimensions = 3 with open(file_name, "r+") as f: #Memory-map the file file_map = mmap.mmap(f.fileno(), 0) position_number = file_map.find(b'NIONS =') file_map.seek(position_number + 7) number_of_atoms = int(file_map.readline()) #Read time step position_number = file_map.find(b'POTIM =') file_map.seek(position_number + 8) time_step = float( file_map.readline().split()[0]) * 1E-3 # in picoseconds #Reading super cell position_number = file_map.find(b'direct lattice vectors') file_map.seek(position_number) file_map.readline() super_cell = [] for i in range(number_of_dimensions): super_cell.append( file_map.readline().split()[0:number_of_dimensions]) super_cell = np.array(super_cell, dtype='double') file_map.seek(position_number) file_map.readline() # Check if number of atoms is multiple of cell atoms if structure is not None: if number_of_atoms % structure.get_number_of_cell_atoms() != 0: print( 'Warning: Number of atoms not matching, check VASP output files' ) # structure.set_number_of_atoms(number_of_atoms) # Read coordinates and energy trajectory = [] energy = [] counter = 0 while True: counter += 1 #Initial cut control if initial_cut > counter: continue position_number = file_map.find(b'POSITION') if position_number < 0: break file_map.seek(position_number) file_map.readline() file_map.readline() read_coordinates = [] for i in range(number_of_atoms): read_coordinates.append( file_map.readline().split()[0:number_of_dimensions]) read_coordinates = np.array(read_coordinates, dtype=float) # in angstrom if template is not None: indexing = np.argsort(template) read_coordinates = read_coordinates[indexing, :] position_number = file_map.find(b'energy(') file_map.seek(position_number) read_energy = file_map.readline().split()[2] trajectory.append(read_coordinates.flatten()) #in angstrom energy.append(np.array(read_energy, dtype=float)) #security routine to limit maximum of steps to read and put in memory if limit_number_steps + initial_cut < counter: print( "Warning! maximum number of steps reached! No more steps will be read" ) break if end_cut is not None and end_cut <= counter: break file_map.close() trajectory = np.array([[[ trajectory[i][j * number_of_dimensions + k] for k in range(number_of_dimensions) ] for j in range(number_of_atoms)] for i in range(len(trajectory))]) if last_steps is not None: trajectory = trajectory[-last_steps:, :, :] energy = energy[-last_steps:] print('Number of total steps read: {0}'.format(trajectory.shape[0])) time = np.array([i * time_step for i in range(trajectory.shape[0])], dtype=float) print('Trajectory file read') return dyn.Dynamics(structure=structure, trajectory=np.array(trajectory, dtype=complex), energy=np.array(energy), time=time, supercell=super_cell, memmap=memmap)
def read_VASP_XDATCAR(file_name, structure=None, time_step=None, limit_number_steps=10000000, last_steps=None, initial_cut=1, end_cut=None, memmap=False, template=None): # Time in picoseconds # Coordinates in Angstroms #Read environtment variables try: temp_directory = os.environ["DYNAPHOPY_TEMPDIR"] if os.path.isdir(temp_directory): print('Set temporal directory: {0}'.format(temp_directory)) temp_directory += '/' else: temp_directory = '' except KeyError: temp_directory = '' number_of_atoms = None bounds = None #Check file exists if not os.path.isfile(file_name): print('Trajectory file does not exist!') exit() #Check time step if time_step is None: print('Warning! XDATCAR file does not contain time step information') print('Using default: 0.001 ps') time_step = 0.001 #Starting reading print("Reading XDATCAR file") print("This could take long, please wait..") #Dimensionality of VASP calculation number_of_dimensions = 3 time = [] data = [] counter = 0 with open(file_name, "r+b") as f: file_map = mmap.mmap(f.fileno(), 0) #Read cell for i in range(2): file_map.readline() a = file_map.readline().split() b = file_map.readline().split() c = file_map.readline().split() super_cell = np.array([a, b, c], dtype='double') for i in range(1): file_map.readline() number_of_atoms = np.array(file_map.readline().split(), dtype=int).sum() while True: counter += 1 #Read time steps position_number = file_map.find(b'Direct configuration') if position_number < 0: break file_map.seek(position_number) time.append(float(file_map.readline().split(b'=')[1])) if memmap: if end_cut: data = np.memmap( temp_directory + 'trajectory.{0}'.format(os.getpid()), dtype='complex', mode='w+', shape=(end_cut - initial_cut + 1, number_of_atoms, number_of_dimensions)) else: print( 'Memory mapping requires to define reading range (use read_from/read_to option)' ) exit() #Initial cut control if initial_cut > counter: continue #Reading coordinates read_coordinates = [] for i in range(number_of_atoms): read_coordinates.append( file_map.readline().split()[0:number_of_dimensions]) read_coordinates = np.array(read_coordinates, dtype=float) # in angstroms if template is not None: indexing = np.argsort(template) read_coordinates = read_coordinates[indexing, :] try: if memmap: data[counter - initial_cut, :, :] = read_coordinates #in angstroms else: data.append(read_coordinates) #in angstroms except ValueError: print("Error reading step {0}".format(counter)) break # print(read_coordinates) #security routine to limit maximum of steps to read and put in memory if limit_number_steps + initial_cut < counter: print( "Warning! maximum number of steps reached! No more steps will be read" ) break if end_cut is not None and end_cut <= counter: break file_map.close() time = np.array(time) * time_step if not memmap: data = np.array(data, dtype=complex) if last_steps is not None: data = data[-last_steps:, :, :] time = time[-last_steps:] return dyn.Dynamics(structure=structure, scaled_trajectory=data, time=time, supercell=super_cell, memmap=memmap)
def read_lammps_trajectory(file_name, structure=None, time_step=None, limit_number_steps=10000000, last_steps=None, initial_cut=1, end_cut=None, memmap=False, template=None): # Time in picoseconds # Coordinates in Angstroms # Read environtment variables try: temp_directory = os.environ["DYNAPHOPY_TEMPDIR"] if os.path.isdir(temp_directory): print('Set temporal directory: {0}'.format(temp_directory)) temp_directory += '/' else: temp_directory = '' except KeyError: temp_directory = '' number_of_atoms = None bounds = None #Check file exists if not os.path.isfile(file_name): print('Trajectory file does not exist!') exit() # Check time step if time_step is None: print( 'Warning! LAMMPS trajectory file does not contain time step information' ) print('Using default: 0.001 ps') time_step = 0.001 # Starting reading print("Reading LAMMPS trajectory") print("This could take long, please wait..") # Dimension of LAMMP calculation if structure is None: number_of_dimensions = 3 else: number_of_dimensions = structure.get_number_of_dimensions() time = [] data = [] counter = 0 lammps_labels = False with open(file_name, "r+") as f: file_map = mmap.mmap(f.fileno(), 0) while True: counter += 1 #Read time steps position_number = file_map.find(b'TIMESTEP') if position_number < 0: break file_map.seek(position_number) file_map.readline() time.append(float(file_map.readline())) if number_of_atoms is None: #Read number of atoms file_map = mmap.mmap(f.fileno(), 0) position_number = file_map.find(b'NUMBER OF ATOMS') file_map.seek(position_number) file_map.readline() number_of_atoms = int(file_map.readline()) # Check if number of atoms is multiple of cell atoms if structure is not None: if number_of_atoms % structure.get_number_of_cell_atoms( ) != 0: print( 'Warning: Number of atoms not matching, check LAMMPS output file' ) if bounds is None: #Read cell file_map = mmap.mmap(f.fileno(), 0) position_number = file_map.find(b'BOX BOUNDS') file_map.seek(position_number) file_map.readline() bounds = [] for i in range(3): bounds.append(file_map.readline().split()) bounds = np.array(bounds, dtype=float) if bounds.shape[1] == 2: bounds = np.append(bounds, np.array([0, 0, 0])[None].T, axis=1) xy = bounds[0, 2] xz = bounds[1, 2] yz = bounds[2, 2] xlo = bounds[0, 0] - np.min([0.0, xy, xz, xy + xz]) xhi = bounds[0, 1] - np.max([0.0, xy, xz, xy + xz]) ylo = bounds[1, 0] - np.min([0.0, yz]) yhi = bounds[1, 1] - np.max([0.0, yz]) zlo = bounds[2, 0] zhi = bounds[2, 1] supercell = np.array([[xhi - xlo, xy, xz], [0, yhi - ylo, yz], [0, 0, zhi - zlo]]).T #for 2D supercell = supercell[:number_of_dimensions, : number_of_dimensions] # Testing cell lx = xhi - xlo ly = yhi - ylo lz = zhi - zlo a = lx b = np.sqrt(pow(ly, 2) + pow(xy, 2)) c = np.sqrt(pow(lz, 2) + pow(xz, 2) + pow(yz, 2)) alpha = np.arccos((xy * xz + ly * yz) / (b * c)) beta = np.arccos(xz / c) gamma = np.arccos(xy / b) # End testing cell # rotate lammps supercell to match unitcell orientation def unit_matrix(matrix): return np.array([ np.array(row) / np.linalg.norm(row) for row in matrix ]) unit_structure = unit_matrix(structure.get_cell()) unit_supercell_lammps = unit_matrix(supercell) transformation_mat = np.dot(np.linalg.inv(unit_structure), unit_supercell_lammps).T supercell = np.dot(supercell, transformation_mat) if memmap: if end_cut: data = np.memmap(temp_directory + 'trajectory.{0}'.format(os.getpid()), dtype='complex', mode='w+', shape=(end_cut - initial_cut + 1, number_of_atoms, number_of_dimensions)) else: print( 'Memory mapping requires to define reading range (use read_from/read_to option)' ) exit() position_number = file_map.find(b'ITEM: ATOMS') file_map.seek(position_number) lammps_labels = file_map.readline() #Initial cut control if initial_cut > counter: time = [] continue #Reading coordinates read_coordinates = [] for i in range(number_of_atoms): read_coordinates.append( file_map.readline().split()[0:number_of_dimensions]) read_coordinates = np.array(read_coordinates, dtype=float) if template is not None: indexing = np.argsort(template) read_coordinates = read_coordinates[indexing, :] try: read_coordinates = np.dot(read_coordinates, transformation_mat) if memmap: data[counter - initial_cut, :, :] = read_coordinates #in angstroms else: data.append(read_coordinates) #in angstroms except ValueError: print("Error reading step {0}".format(counter)) break # print(read_coordinates) #security routine to limit maximum of steps to read and put in memory if limit_number_steps + initial_cut < counter: print( "Warning! maximum number of steps reached! No more steps will be read" ) break if end_cut is not None and end_cut <= counter: break file_map.close() time = np.array(time) * time_step if not memmap: data = np.array(data, dtype=complex) if last_steps is not None: data = data[-last_steps:, :, :] time = time[-last_steps:] # Check position/velocity dump if b'vx vy' in lammps_labels: return dyn.Dynamics(structure=structure, velocity=data, time=time, supercell=supercell, memmap=memmap) if b'x y' in lammps_labels: return dyn.Dynamics(structure=structure, trajectory=data, time=time, supercell=supercell, memmap=memmap) print( 'LAMMPS parsing error. Data not recognized: {}'.format(lammps_labels)) exit()