input/output¶
- class ufedmm.io.StateDataReporter(file, report_interval, **kwargs)[source]¶
Bases:
StateDataReporter
An extension of OpenMM’s app.statedatareporter.StateDataReporter class, which outputs information about a simulation, such as energy, temperature, etc.
All original functionalities are preserved.
Besides, if it is added to an
ExtendedSpaceSimulation
object, e.g. one created through thesimulation()
method, then a new set of keywords are available.- Parameters:
file (str or stream or afed.temperature) – The file to write to, specified as a file name, file object, or
Tee
object.report_interval (int) – The interval (in time steps) at which to report state data.
- Keyword Arguments:
variables (bool, default=False) – Whether to report all dynamical variables related to the extended-space simulation and their associated collective variables.
multipleTemperatures (bool, default=False) – Whether to separately report the temperature estimates for the atoms and for the extended-space dynamical variables.
hillHeights (bool, default=False) – Whether to report the height of the latest deposited metadynamics hill.
collectiveVariables (bool, default=False) – Whether to report the collective variables in all openmm.CustomCVForce objects in the system.
globalParameterStates (pandas.DataFrame, default=None) – A DataFrame containing context global parameters (column names) and sets of values thereof. If it is provided, then the potential energy will be reported for every state defined by these values.
Example
>>> import openmm >>> import ufedmm >>> from openmm import unit >>> from sys import stdout >>> model = ufedmm.AlanineDipeptideModel(water='tip3p') >>> mass = 50*unit.dalton*(unit.nanometer/unit.radians)**2 >>> Ks = 1000*unit.kilojoules_per_mole/unit.radians**2 >>> T = 300*unit.kelvin >>> Ts = 1500*unit.kelvin >>> dt = 2*unit.femtoseconds >>> gamma = 10/unit.picoseconds >>> limit = 180*unit.degrees >>> s_phi = ufedmm.DynamicalVariable( ... 's_phi', -limit, limit, mass, Ts, model.phi, Ks ... ) >>> s_psi = ufedmm.DynamicalVariable( ... 's_psi', -limit, limit, mass, Ts, model.psi, Ks ... ) >>> ufed = ufedmm.UnifiedFreeEnergyDynamics([s_phi, s_psi], T) >>> integrator = ufedmm.GeodesicLangevinIntegrator(T, gamma, dt) >>> platform = openmm.Platform.getPlatformByName('Reference') >>> simulation = ufed.simulation( ... model.topology, model.system, integrator, platform ... ) >>> simulation.context.setPositions(model.positions) >>> simulation.context.setVelocitiesToTemperature(300*unit.kelvin, 1234) >>> reporter = ufedmm.StateDataReporter(stdout, 1, variables=True) >>> reporter.report(simulation, simulation.context.getState()) #"s_phi","phi","s_psi","psi" -3.141592653589793,3.141592653589793,-3.141592653589793,3.141592653589793
- class ufedmm.io.Tee(*files)[source]¶
Bases:
object
Allows the use of multiple outputs in an OpenMM Reporter.
- Parameters:
streams). (A list of valid outputs (file names and/or output) –
Example
>>> import ufedmm >>> import tempfile >>> from sys import stdout >>> file = tempfile.TemporaryFile(mode='w+t') >>> print('test', file=ufedmm.Tee(stdout, file)) test