Skip to content

MM3 Parser

mm3

Parser for Schrödinger MM3* force field files (mm3.fld).

Reads and writes MM3*-format force field files used by Schrödinger's MacroModel, extracting bond, angle, torsion, improper, stretch-bend, and van der Waals parameters for Q2MM optimization.

MM3

MM3(path=None, data=None, method=None, params: list[Param] | None = None, score=None)

Bases: FF

Schrödinger MM3* force field reader and writer (mm3.fld).

Attributes:

Name Type Description
smiles list[str]

MM3* SMILES syntax used in custom parameter sections of the force field file.

sub_names list[str]

Strings used to describe each custom parameter section read.

atom_types list[list[str]]

Atom types derived from each SMILES formula. The SMILES may contain integers referencing earlier atoms, but this list holds only resolved atom type strings.

lines list[str]

Every line from the MM3* force field file.

Initializes an MM3 force field instance.

Parameters:

Name Type Description Default
path str | None

Path to the mm3.fld file. Defaults to None.

None
data object | None

Pre-loaded data. Defaults to None.

None
method str | None

Calculation method label. Defaults to None.

None
params list[Param] | None

Pre-loaded parameters. Defaults to None.

None
score float | None

Objective function score. Defaults to None.

None
Source code in q2mm/parsers/mm3.py
def __init__(self, path=None, data=None, method=None, params: list[Param] | None = None, score=None):
    """Initializes an MM3 force field instance.

    Args:
        path (str | None, optional): Path to the mm3.fld file.
            Defaults to None.
        data (object | None, optional): Pre-loaded data. Defaults to None.
        method (str | None, optional): Calculation method label.
            Defaults to None.
        params (list[Param] | None, optional): Pre-loaded parameters.
            Defaults to None.
        score (float | None, optional): Objective function score.
            Defaults to None.
    """
    super().__init__(path, data, method, params, score)
    self.smiles = []
    self.sub_names = []
    self._atom_types = None
    self._lines = None
    self.atom_type_equivalencies = dict()

atom_types property

atom_types

Derives atom types from SMILES substructure definitions.

Uses the SMILES-esque substructure definition (located directly below the substructure's name) to determine the atom types.

Returns:

Type Description
list[list[str]]

Atom types for each SMILES substructure.

lines property writable

lines

Lines of the force field file.

Returns:

Type Description
list[str]

All lines from the mm3.fld file, lazily loaded from disk on first access.

split_smiles

split_smiles(smiles)

Splits an MM3* SMILES string into individual atom tokens.

Uses the MM3* SMILES substructure definition (located directly below the substructure's name) to determine the atom types.

Parameters:

Name Type Description Default
smiles str

MM3* SMILES substructure string.

required

Returns:

Type Description
list[str]

Individual atom labels extracted from the SMILES string.

Source code in q2mm/parsers/mm3.py
def split_smiles(self, smiles):
    """Splits an MM3* SMILES string into individual atom tokens.

    Uses the MM3* SMILES substructure definition (located directly below the
    substructure's name) to determine the atom types.

    Args:
        smiles (str): MM3* SMILES substructure string.

    Returns:
        (list[str]): Individual atom labels extracted from the SMILES string.
    """
    split_smiles = re.split(co.RE_SPLIT_ATOMS, smiles)
    # I guess this could be an if instead of while since .remove gets rid of
    # all of them, right?
    while "" in split_smiles:
        split_smiles.remove("")
    return split_smiles

convert_smiles_to_types

convert_smiles_to_types(smiles)

Converts an MM3* SMILES string to a list of atom types.

Parameters:

Name Type Description Default
smiles str

MM3* SMILES substructure string.

required

Returns:

Type Description
list[str]

Resolved atom types derived from the SMILES string.

Source code in q2mm/parsers/mm3.py
def convert_smiles_to_types(self, smiles):
    """Converts an MM3* SMILES string to a list of atom types.

    Args:
        smiles (str): MM3* SMILES substructure string.

    Returns:
        (list[str]): Resolved atom types derived from the SMILES string.
    """
    atom_types = self.split_smiles(smiles)
    atom_types = self.convert_to_types(atom_types, atom_types)
    return atom_types

convert_to_types

convert_to_types(atom_labels, atom_types)

Converts atom labels (which may be digit references) to atom types.

For example::

atom_labels = [1, 2]
atom_types  = ["Z0", "P1", "P2"]

would return ["Z0", "P1"].

Parameters:

Name Type Description Default
atom_labels list[str]

Atom labels that can be strings like "C3", "H1", etc. or digit references like "1".

required
atom_types list[str]

Full list of atom type strings to resolve digit references against.

required

Returns:

Type Description
list[str]

Resolved atom types with all digit references replaced by the corresponding atom type string.

Source code in q2mm/parsers/mm3.py
def convert_to_types(self, atom_labels, atom_types):
    """Converts atom labels (which may be digit references) to atom types.

    For example::

        atom_labels = [1, 2]
        atom_types  = ["Z0", "P1", "P2"]

    would return ``["Z0", "P1"]``.

    Args:
        atom_labels (list[str]): Atom labels that can be strings like
            ``"C3"``, ``"H1"``, etc. or digit references like ``"1"``.
        atom_types (list[str]): Full list of atom type strings to
            resolve digit references against.

    Returns:
        (list[str]): Resolved atom types with all digit references
            replaced by the corresponding atom type string.
    """
    return [atom_types[int(x) - 1] if x.strip().isdigit() and x != "00" else x for x in atom_labels]

import_ff

import_ff(path=None, sub_search='OPT')

Reads parameters from an mm3.fld file.

Parses bond, angle, stretch-bend, torsion, improper, and van der Waals parameters from the MM3* force field file.

Parameters:

Name Type Description Default
path str | None

Path to the force field file. Defaults to self.path.

None
sub_search str

Marker string used to identify optimizable substructure sections. Defaults to "OPT".

'OPT'
Source code in q2mm/parsers/mm3.py
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
def import_ff(self, path=None, sub_search="OPT"):
    """Reads parameters from an mm3.fld file.

    Parses bond, angle, stretch-bend, torsion, improper, and van der
    Waals parameters from the MM3* force field file.

    Args:
        path (str | None, optional): Path to the force field file.
            Defaults to ``self.path``.
        sub_search (str, optional): Marker string used to identify
            optimizable substructure sections. Defaults to ``"OPT"``.
    """
    if path is None:
        path = self.path
    self.params: list[Param] = []
    self.smiles = []
    self.sub_names = []
    with open(path) as f:
        logger.log(15, f"READING: {path}")
        section_sub = False
        section_smiles = False
        section_vdw = False
        section_atm_eqv = False  # atom type equivalencies section
        for i, line in enumerate(f):
            if section_atm_eqv:
                if line.startswith(" C") and len(self.atom_type_equivalencies.items()) > 0:
                    section_atm_eqv = False
                    continue
                elif not line.startswith(" C") and not line.startswith("-5"):
                    equivalency = [typ.strip() for typ in line.split()[1:]]
                    for typ in equivalency[1:]:
                        self.atom_type_equivalencies[typ] = equivalency[0]
                    continue

            # These lines are for parameters.
            if not section_sub and sub_search in line and line.startswith(" C"):
                matched = re.match(rf"\sC\s+({co.RE_SUB})\s+", line)
                assert matched is not None, f"[L{i + 1}] Can't read substructure name: {line}"
                if matched is not None:
                    # Oh good, you found your substructure!
                    section_sub = True
                    sub_name = matched.group(1).strip()
                    self.sub_names.append(sub_name)
                    logger.log(
                        15,
                        f"[L{i + 1}] Start of substructure: {sub_name}",
                    )
                    section_smiles = True
                    continue
            elif section_smiles is True:
                matched = re.match(rf"\s9\s+({co.RE_SMILES})\s", line)
                assert matched is not None, f"[L{i + 1}] Can't read substructure SMILES: {line}"
                smiles = matched.group(1)
                self.smiles.append(smiles)
                logger.log(15, f"  -- SMILES: {self.smiles[-1]}")
                logger.log(15, "  -- Atom types: {}".format(" ".join(self.atom_types[-1])))
                section_smiles = False
                continue
            # Marks the end of a substructure.
            elif section_sub and line.startswith("-3"):
                logger.log(
                    15,
                    f"[L{i}] End of substructure: {self.sub_names[-1]}",
                )
                section_sub = False
                continue
            if "OPT" in line and section_vdw:
                logger.log(
                    5,
                    "[L{}] Found Van der Waals:\n{}".format(i + 1, line.strip("\n")),
                )
                atm = line[2:5]
                rad = line[5:15]
                eps = line[16:26]
                self.params.extend(
                    (
                        Param(
                            atom_types=atm,
                            ptype="vdwr",
                            ff_col=1,
                            ff_row=i + 1,
                            value=float(rad),
                        ),
                        Param(
                            atom_types=atm,
                            ptype="vdwe",
                            ff_col=2,
                            ff_row=i + 1,
                            value=float(eps),
                        ),
                    )
                )
                continue
            if "OPT" in line or section_sub:
                # Bonds.
                if match_mm3_bond(line):
                    logger.log(5, "[L{}] Found bond:\n{}".format(i + 1, line.strip("\n")))
                    if section_sub:
                        atm_lbls = [line[4:6], line[8:10]]
                        atm_typs = self.convert_to_types(atm_lbls, self.atom_types[-1])
                    else:
                        atm_typs = [line[4:6], line[9:11]]
                        atm_lbls = atm_typs
                        comment = line[COM_POS_START:].strip()
                        self.sub_names.append(comment)
                    parm_cols = line[P_1_START:P_3_END]
                    parm_cols = [float(x) for x in parm_cols.split()]
                    self.params.extend(
                        (
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="be",
                                ff_col=1,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[0],
                            ),
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="bf",
                                ff_col=2,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[1],
                            ),
                        )
                    )
                    # Some bonds parameters don't use bond dipoles.
                    with contextlib.suppress(IndexError):
                        self.params.append(
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="q",
                                ff_col=3,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[2],
                            )
                        )
                    continue
                # Angles.
                elif match_mm3_angle(line):
                    logger.log(5, "[L{}] Found angle:\n{}".format(i + 1, line.strip("\n")))
                    if section_sub:
                        # Do stuff.
                        atm_lbls = [line[4:6], line[8:10], line[12:14]]
                        atm_typs = self.convert_to_types(atm_lbls, self.atom_types[-1])
                    else:
                        # Do other method.
                        atm_typs = [line[4:6], line[9:11], line[14:16]]
                        atm_lbls = atm_typs
                        comment = line[COM_POS_START:].strip()
                        self.sub_names.append(comment)
                    parm_cols = line[P_1_START:P_3_END]
                    parm_cols = [float(x) for x in parm_cols.split()]
                    self.params.extend(
                        (
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="ae",
                                ff_col=1,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[0],
                            ),
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="af",
                                ff_col=2,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[1],
                            ),
                        )
                    )
                    continue
                # Stretch-bends.
                elif match_mm3_stretch_bend(line):
                    logger.log(
                        5,
                        "[L{}] Found stretch-bend:\n{}".format(i + 1, line.strip("\n")),
                    )
                    if section_sub:
                        # Do stuff.
                        atm_lbls = [line[4:6], line[8:10], line[12:14]]
                        atm_typs = self.convert_to_types(atm_lbls, self.atom_types[-1])
                    else:
                        # Do other method.
                        atm_typs = [line[4:6], line[9:11], line[14:16]]
                        atm_lbls = atm_typs
                        comment = line[COM_POS_START:].strip()
                        self.sub_names.append(comment)
                    parm_cols = line[P_1_START:P_3_END]
                    parm_cols = [float(x) for x in parm_cols.split()]
                    self.params.append(
                        Param(
                            atom_labels=atm_lbls,
                            atom_types=atm_typs,
                            ptype="sb",
                            ff_col=1,
                            ff_row=i + 1,
                            label=line[:2],
                            value=parm_cols[0],
                        )
                    )
                    continue
                # Torsions.
                elif match_mm3_lower_torsion(line):
                    logger.log(
                        5,
                        "[L{}] Found torsion:\n{}".format(i + 1, line.strip("\n")),
                    )
                    if section_sub:
                        # Do stuff.
                        atm_lbls = [line[4:6], line[8:10], line[12:14], line[16:18]]
                        atm_typs = self.convert_to_types(atm_lbls, self.atom_types[-1])
                    else:
                        # Do other method.
                        atm_typs = [line[4:6], line[9:11], line[14:16], line[19:21]]
                        atm_lbls = atm_typs
                        comment = line[COM_POS_START:].strip()
                        self.sub_names.append(comment)
                    parm_cols = line[P_1_START:P_3_END]
                    parm_cols = [float(x) for x in parm_cols.split()]
                    self.params.extend(
                        (
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="df",
                                ff_col=1,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[0],
                            ),
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="df",
                                ff_col=2,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[1],
                            ),
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="df",
                                ff_col=3,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[2],
                            ),
                        )
                    )
                    continue
                # Higher order torsions.
                elif match_mm3_higher_torsion(line):
                    logger.log(
                        5,
                        "[L{}] Found higher order torsion:\n{}".format(i + 1, line.strip("\n")),
                    )
                    # Will break if torsions aren't also looked up.
                    atm_lbls = self.params[-1].atom_labels
                    atm_typs = self.params[-1].atom_types
                    parm_cols = line[P_1_START:P_3_END]
                    parm_cols = [float(x) for x in parm_cols.split()]
                    self.params.extend(
                        (
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="df",
                                ff_col=1,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[0],
                            ),
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="df",
                                ff_col=2,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[1],
                            ),
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="df",
                                ff_col=3,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[2],
                            ),
                        )
                    )
                    continue
                # Improper torsions.
                elif match_mm3_improper(line):
                    logger.log(
                        5,
                        "[L{}] Found torsion:\n{}".format(i + 1, line.strip("\n")),
                    )
                    if section_sub:
                        # Do stuff.
                        atm_lbls = [line[4:6], line[8:10], line[12:14], line[16:18]]
                        atm_typs = self.convert_to_types(atm_lbls, self.atom_types[-1])
                    else:
                        # Do other method.
                        atm_typs = [line[4:6], line[9:11], line[14:16], line[19:21]]
                        atm_lbls = atm_typs
                        comment = line[COM_POS_START:].strip()
                        self.sub_names.append(comment)
                    parm_cols = line[P_1_START:P_3_END]
                    parm_cols = [float(x) for x in parm_cols.split()]
                    self.params.extend(
                        (
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="imp1",
                                ff_col=1,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[0],
                            ),
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="imp2",
                                ff_col=2,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[1],
                            ),
                        )
                    )
                    continue
                # Bonds.
                elif match_mm3_vdw(line):
                    logger.log(5, "[L{}] Found vdw:\n{}".format(i + 1, line.strip("\n")))
                    if section_sub:
                        atm_lbls = [line[4:6], line[8:10]]
                        atm_typs = self.convert_to_types(atm_lbls, self.atom_types[-1])
                    parm_cols = line[P_1_START:P_3_END]
                    parm_cols = [float(x) for x in parm_cols.split()]
                    self.params.extend(
                        (
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="vdwr",
                                ff_col=1,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[0],
                            ),
                            Param(
                                atom_labels=atm_lbls,
                                atom_types=atm_typs,
                                ptype="vdwfc",
                                ff_col=2,
                                ff_row=i + 1,
                                label=line[:2],
                                value=parm_cols[1],
                            ),
                        )
                    )
                    continue
            # The Van der Waals are stored in annoying way.
            if line.startswith("-6"):
                section_vdw = True
                continue
            if "New Atom Type Equivalencies" in line:
                section_atm_eqv = True
                continue
    logger.log(15, f"  -- Read {len(self.params)} parameters.")

export_ff

export_ff(path=None, params=None, lines=None)

Exports the force field to a file, typically mm3.fld.

Parameters:

Name Type Description Default
path str | None

Output file path. Defaults to self.path.

None
params list[Param] | None

Parameters to write. Defaults to self.params.

None
lines list[str] | None

Base file lines to modify. Generated via readlines() on the mm3.fld file. Defaults to self.lines.

None
Source code in q2mm/parsers/mm3.py
def export_ff(self, path=None, params=None, lines=None):
    """Exports the force field to a file, typically mm3.fld.

    Args:
        path (str | None, optional): Output file path. Defaults to
            ``self.path``.
        params (list[Param] | None, optional): Parameters to write.
            Defaults to ``self.params``.
        lines (list[str] | None, optional): Base file lines to modify.
            Generated via ``readlines()`` on the mm3.fld file.
            Defaults to ``self.lines``.
    """
    if path is None:
        path = self.path
    if params is None:
        params: list[Param] = self.params
    if lines is None:
        lines = self.lines
    for param in params:
        logger.log(1, f">>> param: {param} param.value: {param.value}")
        line = lines[param.ff_row - 1]
        # There are some problems with this. Probably an optimization
        # technique gave you these crazy parameter values. Ideally, this
        # entire trial FF should be discarded.
        # Someday export_ff should raise an exception when these values
        # get too rediculous, and this exception should be handled by the
        # optimization techniques appropriately.
        if abs(param.value) > 999.0:
            logger.warning(f"Value of {param} is too high! Skipping write.")
        elif param.ff_col == 1:
            lines[param.ff_row - 1] = line[:P_1_START] + f"{param.value:10.4f}" + line[P_1_END:]
        elif param.ff_col == 2:
            lines[param.ff_row - 1] = line[:P_2_START] + f"{param.value:10.4f}" + line[P_2_END:]
        elif param.ff_col == 3:
            lines[param.ff_row - 1] = line[:P_3_START] + f"{param.value:10.4f}" + line[P_3_END:]
    with open(path, "w") as f:
        f.writelines(lines)
    logger.log(10, f"WROTE: {path}")

get_DOFs_by_ff_row

get_DOFs_by_ff_row(structs: list[Structure]) -> dict

Groups degrees of freedom by force field row number.

Parameters:

Name Type Description Default
structs list[Structure]

Structures whose DOFs are collected.

required

Returns:

Type Description
dict

Mapping of ff_row to a list of :class:DOF instances (bonds, angles, and torsions).

Source code in q2mm/parsers/mm3.py
def get_DOFs_by_ff_row(self, structs: list[Structure]) -> dict:
    """Groups degrees of freedom by force field row number.

    Args:
        structs (list[Structure]): Structures whose DOFs are collected.

    Returns:
        (dict): Mapping of ``ff_row`` to a list of :class:`DOF` instances
            (bonds, angles, and torsions).
    """
    dof_by_param = dict()
    for param in self.params:
        dof_by_param[param.ff_row]: list[DOF] = []
    for struct in structs:
        for bond in struct.bonds:
            dof_by_param[bond.ff_row].append(bond)
        for angle in struct.angles:
            dof_by_param[angle.ff_row].append(angle)
        for dihed in struct.torsions:
            dof_by_param[dihed.ff_row].append(dihed)
    return dof_by_param

match_mm3_label

match_mm3_label(mm3_label)

Checks whether a line has a recognized MM3* parameter label.

The label is the first 2 characters in the line containing the parameter in a Schrödinger mm3.fld file.

Parameters:

Name Type Description Default
mm3_label str

Line or string whose first 2 characters are checked.

required

Returns:

Type Description
Match | None

Match object if the label is recognized, else None.

Source code in q2mm/parsers/mm3.py
def match_mm3_label(mm3_label):
    """Checks whether a line has a recognized MM3* parameter label.

    The label is the first 2 characters in the line containing the parameter
    in a Schrödinger mm3.fld file.

    Args:
        mm3_label (str): Line or string whose first 2 characters are checked.

    Returns:
        (re.Match | None): Match object if the label is recognized, else None.
    """
    return re.match(r"[\s5a-z][1-5]", mm3_label)

match_mm3_vdw

match_mm3_vdw(mm3_label)

Matches MM3* label for van der Waals parameters.

Parameters:

Name Type Description Default
mm3_label str

Line or string whose first 2 characters are checked.

required

Returns:

Type Description
Match | None

Match object if the label matches, else None.

Source code in q2mm/parsers/mm3.py
def match_mm3_vdw(mm3_label):
    """Matches MM3* label for van der Waals parameters.

    Args:
        mm3_label (str): Line or string whose first 2 characters are checked.

    Returns:
        (re.Match | None): Match object if the label matches, else None.
    """
    return re.match(r"[\sa-z]6", mm3_label)

match_mm3_bond

match_mm3_bond(mm3_label)

Matches MM3* label for bonds.

Parameters:

Name Type Description Default
mm3_label str

Line or string whose first 2 characters are checked.

required

Returns:

Type Description
Match | None

Match object if the label matches, else None.

Source code in q2mm/parsers/mm3.py
def match_mm3_bond(mm3_label):
    """Matches MM3* label for bonds.

    Args:
        mm3_label (str): Line or string whose first 2 characters are checked.

    Returns:
        (re.Match | None): Match object if the label matches, else None.
    """
    return re.match(r"[\sa-z]1", mm3_label)

match_mm3_angle

match_mm3_angle(mm3_label)

Matches MM3* label for angles.

Parameters:

Name Type Description Default
mm3_label str

Line or string whose first 2 characters are checked.

required

Returns:

Type Description
Match | None

Match object if the label matches, else None.

Source code in q2mm/parsers/mm3.py
def match_mm3_angle(mm3_label):
    """Matches MM3* label for angles.

    Args:
        mm3_label (str): Line or string whose first 2 characters are checked.

    Returns:
        (re.Match | None): Match object if the label matches, else None.
    """
    return re.match(r"[\sa-z]2", mm3_label)

match_mm3_stretch_bend

match_mm3_stretch_bend(mm3_label)

Matches MM3* label for stretch-bends.

Parameters:

Name Type Description Default
mm3_label str

Line or string whose first 2 characters are checked.

required

Returns:

Type Description
Match | None

Match object if the label matches, else None.

Source code in q2mm/parsers/mm3.py
def match_mm3_stretch_bend(mm3_label):
    """Matches MM3* label for stretch-bends.

    Args:
        mm3_label (str): Line or string whose first 2 characters are checked.

    Returns:
        (re.Match | None): Match object if the label matches, else None.
    """
    return re.match(r"[\sa-z]3", mm3_label)

match_mm3_torsion

match_mm3_torsion(mm3_label)

Matches MM3* label for all orders of torsional parameters.

Parameters:

Name Type Description Default
mm3_label str

Line or string whose first 2 characters are checked.

required

Returns:

Type Description
Match | None

Match object if the label matches, else None.

Source code in q2mm/parsers/mm3.py
def match_mm3_torsion(mm3_label):
    """Matches MM3* label for all orders of torsional parameters.

    Args:
        mm3_label (str): Line or string whose first 2 characters are checked.

    Returns:
        (re.Match | None): Match object if the label matches, else None.
    """
    return re.match(r"[\sa-z]4|54", mm3_label)

match_mm3_lower_torsion

match_mm3_lower_torsion(mm3_label)

Matches MM3* label for torsions (1st through 3rd order).

Parameters:

Name Type Description Default
mm3_label str

Line or string whose first 2 characters are checked.

required

Returns:

Type Description
Match | None

Match object if the label matches, else None.

Source code in q2mm/parsers/mm3.py
def match_mm3_lower_torsion(mm3_label):
    """Matches MM3* label for torsions (1st through 3rd order).

    Args:
        mm3_label (str): Line or string whose first 2 characters are checked.

    Returns:
        (re.Match | None): Match object if the label matches, else None.
    """
    return re.match(r"[\sa-z]4", mm3_label)

match_mm3_higher_torsion

match_mm3_higher_torsion(mm3_label)

Matches MM3* label for torsions (4th through 6th order).

Parameters:

Name Type Description Default
mm3_label str

Line or string whose first 2 characters are checked.

required

Returns:

Type Description
Match | None

Match object if the label matches, else None.

Source code in q2mm/parsers/mm3.py
def match_mm3_higher_torsion(mm3_label):
    """Matches MM3* label for torsions (4th through 6th order).

    Args:
        mm3_label (str): Line or string whose first 2 characters are checked.

    Returns:
        (re.Match | None): Match object if the label matches, else None.
    """
    return re.match("54", mm3_label)

match_mm3_improper

match_mm3_improper(mm3_label)

Matches MM3* label for improper torsions.

Parameters:

Name Type Description Default
mm3_label str

Line or string whose first 2 characters are checked.

required

Returns:

Type Description
Match | None

Match object if the label matches, else None.

Source code in q2mm/parsers/mm3.py
def match_mm3_improper(mm3_label):
    """Matches MM3* label for improper torsions.

    Args:
        mm3_label (str): Line or string whose first 2 characters are checked.

    Returns:
        (re.Match | None): Match object if the label matches, else None.
    """
    return re.match(r"[\sa-z]5", mm3_label)