@@ -190,7 +190,8 @@ def compute_vaf(
190190 fft : bool = False ,
191191 index : SliceLike = (0 , None , 1 ),
192192 filter_atoms : MaybeSequence [MaybeSequence [Optional [int ]]] = ((None ),),
193- ) -> NDArray [float64 ]:
193+ time_step : float = 1.0 ,
194+ ) -> tuple [NDArray [float64 ], list [NDArray [float64 ]]]:
194195 """
195196 Compute the velocity autocorrelation function (VAF) of `data`.
196197
@@ -212,11 +213,36 @@ def compute_vaf(
212213 filter_atoms : MaybeSequence[MaybeSequence[Optional[int]]]
213214 Compute the VAF averaged over subsets of the system.
214215 Default is all atoms.
216+ time_step : float
217+ Time step for scaling lags to align with input data.
218+ Default is 1 (i.e. no scaling).
215219
216220 Returns
217221 -------
218- MaybeSequence[NDArray[float64]]
222+ lags : numpy.ndarray
223+ Lags at which the VAFs have been computed.
224+ vafs : list[numpy.ndarray]
219225 Computed VAF(s).
226+
227+ Notes
228+ -----
229+ `filter_atoms` is given as a series of sequences of atoms, where
230+ each element in the series denotes a VAF subset to calculate and
231+ each sequence determines the atoms (by index) to be included in that VAF.
232+
233+ E.g.
234+
235+ .. code-block: Python
236+
237+ # Species indices in cell
238+ na = (1, 3, 5, 7)
239+ cl = (2, 4, 6, 8)
240+
241+ compute_vaf(..., filter_atoms=(na, cl))
242+
243+ Would compute separate VAFs for each species.
244+
245+ By default, one VAF will be computed for all atoms in the structure.
220246 """
221247 # Ensure if passed scalars they are turned into correct dimensionality
222248 if not isinstance (filter_atoms , Sequence ):
@@ -270,17 +296,24 @@ def compute_vaf(
270296
271297 vafs /= n_steps - np .arange (n_steps )
272298
299+ lags = np .arange (n_steps ) * time_step
300+
273301 if fft :
274302 vafs = np .fft .fft (vafs , axis = 0 )
275-
276- vafs = [
277- np .average ([vafs [used_atoms [i ]] for i in atoms ], axis = 0 )
278- for atoms in filter_atoms
279- ]
303+ lags = np .fft .fftfreq (n_steps , time_step )
304+
305+ vafs = (
306+ lags ,
307+ [
308+ np .average ([vafs [used_atoms [i ]] for i in atoms ], axis = 0 )
309+ for atoms in filter_atoms
310+ ],
311+ )
280312
281313 if filenames :
282- for filename , vaf in zip (filenames , vafs ):
314+ for vaf , filename in zip (vafs [ 1 ], filenames ):
283315 with open (filename , "w" , encoding = "utf-8" ) as out_file :
284- print (* vaf , file = out_file , sep = "\n " )
316+ for lag , dat in zip (lags , vaf ):
317+ print (lag , dat , file = out_file )
285318
286319 return vafs
0 commit comments