def particles(self, filename='pytim_part.vtk', group=None, sequence=False): """ Write to vtk files the particles in a group. :param str filename: the file name :param bool sequence: if true writes a sequence of files adding the frame to the filename :param AtomGroup group: if None, writes the whole universe >>> import MDAnalysis as mda >>> import pytim >>> from pytim.datafiles import MICELLE_PDB >>> u = mda.Universe(MICELLE_PDB) >>> g = u.select_atoms('resname DPC') >>> inter= pytim.WillardChandler(u, group=g, alpha=3.0, mesh=2.0) >>> # writes on part.vtk >>> inter.writevtk.particles('part.vtk') >>> # writes on part.<frame>.vtk >>> inter.writevtk.particles('part.vtk',sequence=True) """ inter = self.interface if sequence is True: filename = vtk.consecutive_filename(inter.universe, filename) if group is None: group = inter.universe.atoms self._dump_group(group, filename)
def density(self, filename='pytim_dens.vtk', sequence=False): """ Write to vtk files the volumetric density: :param str filename: the file name :param bool sequence: if true writes a sequence of files adding the frame to the filename >>> import MDAnalysis as mda >>> import pytim >>> from pytim.datafiles import MICELLE_PDB >>> u = mda.Universe(MICELLE_PDB) >>> g = u.select_atoms('resname DPC') >>> interface = pytim.WillardChandler(u, group=g, alpha=3.0, mesh=2.0) >>> interface.writevtk.density('dens.vtk') # writes on dens.vtk >>> interface.writevtk.density('dens.vtk',sequence=True) # writes on dens.<frame>.vtk """ inter = self.interface if sequence == True: filename = vtk.consecutive_filename(inter.universe, filename) vtk.write_scalar_grid(filename, inter.ngrid, inter.spacing, inter.density_field)
def surface(self, filename='pytim_surf.vtk', sequence=False): """ Write to vtk files the triangulated surface: :param str filename: the file name :param bool sequence: if true writes a sequence of files adding the frame to the filename >>> import MDAnalysis as mda >>> import pytim >>> from pytim.datafiles import MICELLE_PDB >>> u = mda.Universe(MICELLE_PDB) >>> g = u.select_atoms('resname DPC') >>> interface = pytim.WillardChandler(u, group=g, alpha=3.0, mesh=2.0) >>> interface.writevtk.surface('surf.vtk') # writes on surf.vtk >>> interface.writevtk.surface('surf.vtk',sequence=True) # writes on surf.<frame>.vtk """ inter = self.interface vertices = inter.triangulated_surface[0] faces = inter.triangulated_surface[1] normals = inter.triangulated_surface[2] if sequence == True: filename = vtk.consecutive_filename(inter.universe, filename) vtk.write_triangulation(filename, vertices[::, ::-1], faces, normals)