Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

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

538

539

540

541

542

543

544

545

546

547

548

549

550

551

552

553

554

555

556

557

558

559

560

561

562

563

564

565

566

567

568

569

570

571

572

573

574

575

576

577

578

579

580

581

582

583

584

585

586

587

588

589

590

591

592

593

594

595

596

597

598

599

600

601

602

603

604

605

606

607

608

609

610

611

612

613

614

615

616

617

618

619

620

621

622

623

624

625

626

627

628

629

630

631

632

633

634

635

636

637

638

639

640

641

642

643

644

645

646

647

648

649

650

651

652

653

654

655

656

657

658

659

660

661

662

663

664

665

666

667

668

669

670

671

672

673

674

675

676

677

678

679

680

681

682

683

684

685

686

687

688

689

690

691

692

693

694

695

696

697

698

699

700

701

702

703

704

705

706

707

708

709

710

711

712

713

714

715

716

717

718

719

720

721

722

723

724

725

726

727

728

729

730

731

732

733

734

735

736

737

738

739

740

# coding: utf-8 

# Copyright (c) Pymatgen Development Team. 

# Distributed under the terms of the MIT License. 

 

from __future__ import division, print_function, unicode_literals, \ 

absolute_import 

 

""" 

This module implements classes for generating/parsing Lammps data files. 

""" 

 

from six.moves import range 

from io import open 

import re 

from collections import OrderedDict 

 

from monty.json import MSONable, MontyDecoder 

 

from pymatgen.core.structure import Molecule, Structure 

 

__author__ = 'Kiran Mathew' 

__email__ = "kmathew@lbl.gov" 

 

 

class LammpsData(MSONable): 

""" 

Basic Lammps data: just the atoms section 

 

Args: 

box_size (list): [[x_min, x_max], [y_min,y_max], [z_min,z_max]] 

atomic_masses (list): [[atom type, mass],...] 

atoms_data (list): 

[[atom id, mol id, atom type, charge, x, y, z ...], ... ] 

""" 

 

def __init__(self, box_size, atomic_masses, atoms_data): 

self.box_size = box_size 

self.natoms = len(atoms_data) 

self.natom_types = len(atomic_masses) 

self.atomic_masses = list(atomic_masses) 

self.atoms_data = atoms_data 

 

def __str__(self): 

""" 

string representation of LammpsData 

 

Returns: 

String representation of the data file 

""" 

lines = [] 

lines.append("Data file generated by rubicon\n") 

lines.append("{} atoms\n".format(self.natoms)) 

lines.append("{} atom types\n".format(self.natom_types)) 

lines.append("{} {} xlo xhi\n{} {} ylo yhi\n{} {} zlo zhi".format( 

self.box_size[0][0], self.box_size[0][1], 

self.box_size[1][0], self.box_size[1][1], 

self.box_size[2][0], self.box_size[2][1])) 

self.set_lines_from_list(lines, "Masses", self.atomic_masses) 

self.set_lines_from_list(lines, "Atoms", self.atoms_data) 

return '\n'.join(lines) 

 

def write_data_file(self, filename): 

""" 

write lammps data input file from the string representation 

of the data. 

 

Args: 

filename (string): data file name 

""" 

with open(filename, 'w') as f: 

f.write(self.__str__()) 

 

@staticmethod 

def get_basic_system_info(structure): 

""" 

Return basic system info from the given structure 

 

Args: 

structure (Structure) 

 

Returns: 

number of atoms, number of atom types, box size, mapping 

between the atom id and corresponding atomic masses 

""" 

natoms = len(structure) 

natom_types = len(structure.symbol_set) 

box_size = [[0.0, structure.lattice.a], 

[0.0, structure.lattice.b], 

[0.0, structure.lattice.c]] 

elements = structure.composition.elements 

elements = sorted(elements, key=lambda el: el.atomic_mass) 

atomic_masses_dict = OrderedDict( 

[(el.symbol, [i + 1, el.data["Atomic mass"]]) 

for i, el in enumerate(elements)]) 

return natoms, natom_types, box_size, atomic_masses_dict 

 

@staticmethod 

def get_atoms_data(structure, atomic_masses_dict, set_charge=True): 

""" 

return the atoms data: 

atom_id, molecule tag, atom_type, charge(if present else 0), x, y, z. 

The molecule_tag is set to 1(i.e the whole structure corresponds to 

just one molecule). This corresponds to lammps command: "atom_style 

charge" 

 

Args: 

structure (Structure) 

atomic_masses_dict (dict): 

{ atom symbol : [atom_id, atomic mass], ... } 

set_charge (bool): whether or not to set the charge field in Atoms 

 

Returns: 

[[atom_id, molecule tag, atom_type, charge(if present), x, y, z], 

... ] 

""" 

atoms_data = [] 

for i, site in enumerate(structure): 

atom_type = atomic_masses_dict[site.specie.symbol][0] 

if set_charge: 

if hasattr(site, "charge"): 

atoms_data.append([i + 1, 1, atom_type, site.charge, 

site.x, site.y, site.z]) 

else: 

atoms_data.append([i + 1, 1, atom_type, 0.0, 

site.x, site.y, site.z]) 

else: 

atoms_data.append( 

[i + 1, 1, atom_type, site.x, site.y, site.z]) 

return atoms_data 

 

@staticmethod 

def set_lines_from_list(lines, block_name, input_list): 

""" 

Append the values from the input list that corresponds to the block 

with name 'block_name' to the list of lines. 

 

Args: 

lines (list) 

block_name (string): name of the data block, e.g. 'Atoms', 

'Bonds' etc 

input_list (list): list of values 

""" 

if input_list: 

lines.append("\n" + block_name + " \n") 

_ = [lines.append(" ".join([str(x) for x in ad])) for ad in 

input_list] 

 

@staticmethod 

def from_structure(input_structure, box_size, set_charge=True): 

""" 

Set LammpsData from the given structure. If the input structure is 

a Molecule, the box_size will be used to created a boxed molecule 

else the sites within the structure will be set in a box with 

dimension set from either box_size or the lattice of the structure, 

whichever is bigger. 

 

Args: 

input_structure (Structure) 

box_size (list): [[x_min, x_max], [y_min, y_max], [z_min, z_max]] 

set_charge (bool): whether or not to set the charge field in 

Atoms. If true, the charge will be non-zero only if the 

input_structure has the "charge" site property set. 

 

Returns: 

LammpsData 

""" 

boxed_molecule = get_boxed_molecule(input_structure, box_size) 

natoms, natom_types, box_size, atomic_masses_dict = \ 

LammpsData.get_basic_system_info(boxed_molecule) 

atoms_data = LammpsData.get_atoms_data(boxed_molecule, 

atomic_masses_dict, 

set_charge=set_charge) 

return LammpsData(box_size, atomic_masses_dict.values(), atoms_data) 

 

@staticmethod 

def from_file(data_file, read_charge=True): 

""" 

Return LammpsData object from the data file. 

Note: use this to read in data files that conform with 

atom_style = charge or atomic 

 

Args: 

data_file (string): data file name 

read_charge (bool): if true, read in data files that conform with 

atom_style = charge else atom_style = atomic 

 

Returns: 

LammpsData 

""" 

atomic_masses = [] # atom_type(starts from 1): mass 

box_size = [] 

atoms_data = [] 

# atom_id, mol_id, atom_type, charge, x, y, z 

if read_charge: 

atoms_pattern = re.compile( 

"^\s*(\d+)\s+(\d+)\s+(\d+)\s+([0-9eE\.+-]+)\s+(" 

"[0-9eE\.+-]+)\s+([0-9eE\.+-]+)\s+([0-9eE\.+-]+)\w*") 

# atom_id, mol_id, atom_type, x, y, z 

else: 

atoms_pattern = re.compile( 

"^\s*(\d+)\s+(\d+)\s+(\d+)\s+([0-9eE\.+-]+)\s+(" 

"[0-9eE\.+-]+)\s+([0-9eE\.+-]+)\w*") 

# atom_type, mass 

masses_pattern = re.compile("^\s*(\d+)\s+([0-9\.]+)$") 

box_pattern = re.compile( 

"^([0-9eE\.+-]+)\s+([0-9eE\.+-]+)\s+[xyz]lo\s+[xyz]hi") 

with open(data_file) as df: 

for line in df: 

if masses_pattern.search(line): 

m = masses_pattern.search(line) 

atomic_masses.append([int(m.group(1)), float(m.group(2))]) 

if box_pattern.search(line): 

m = box_pattern.search(line) 

box_size.append([float(m.group(1)), float(m.group(2))]) 

m = atoms_pattern.search(line) 

if m: 

# atom id, mol id, atom type 

line_data = [int(i) for i in m.groups()[:3]] 

# charge, x, y, z 

line_data.extend([float(i) for i in m.groups()[3:]]) 

atoms_data.append(line_data) 

return LammpsData(box_size, atomic_masses, atoms_data) 

 

def as_dict(self): 

d = MSONable.as_dict(self) 

if hasattr(self, "kwargs"): 

d.update(**self.kwargs) 

return d 

 

@classmethod 

def from_dict(cls, d): 

decoded = {k: MontyDecoder().process_decoded(v) for k, v in d.items() 

if not k.startswith("@")} 

return cls(**decoded) 

 

 

class LammpsForceFieldData(LammpsData): 

""" 

Sets Lammps data input file from force field parameters. 

 

Args: 

box_size (list): [[x_min,x_max], [y_min,y_max], [z_min,z_max]] 

atomic_masses (list): [ [atom type, atomic mass], ... ] 

pair_coeffs (list): vdw coefficients, 

[[unique id, sigma, epsilon ], ... ] 

bond_coeffs (list): bond coefficients, 

[[unique id, value1, value2 ], ... ] 

angle_coeffs (list): angle coefficients, 

[[unique id, value1, value2, value3 ], ... ] 

dihedral_coeffs (list): dihedral coefficients, 

[[unique id, value1, value2, value3, value4], ... ] 

improper_coeffs (list): improper dihedral coefficients, 

[[unique id, value1, value2, value3, value4], ... ] 

atoms_data (list): [[atom id, mol id, atom type, charge, x,y,z, ...], ... ] 

bonds_data (list): [[bond id, bond type, value1, value2], ... ] 

angles_data (list): [[angle id, angle type, value1, value2, value3], ... ] 

dihedrals_data (list): 

[[dihedral id, dihedral type, value1, value2, value3, value4], ... ] 

imdihedrals_data (list): 

[[improper dihedral id, improper dihedral type, value1, value2, 

value3, value4], ... ] 

""" 

 

def __init__(self, box_size, atomic_masses, pair_coeffs, bond_coeffs, 

angle_coeffs, dihedral_coeffs, improper_coeffs, atoms_data, 

bonds_data, angles_data, dihedrals_data, imdihedrals_data): 

super(LammpsForceFieldData, self).__init__(box_size, atomic_masses, 

atoms_data) 

# number of types 

self.nbond_types = len(bond_coeffs) 

self.nangle_types = len(angle_coeffs) 

self.ndih_types = len(dihedral_coeffs) 

self.nimdih_types = len(improper_coeffs) 

# number of parameters 

self.nbonds = len(bonds_data) 

self.nangles = len(angles_data) 

self.ndih = len(dihedrals_data) 

self.nimdihs = len(imdihedrals_data) 

# coefficients 

self.pair_coeffs = pair_coeffs 

self.bond_coeffs = bond_coeffs 

self.angle_coeffs = angle_coeffs 

self.dihedral_coeffs = dihedral_coeffs 

self.improper_coeffs = improper_coeffs 

# data 

self.bonds_data = bonds_data 

self.angles_data = angles_data 

self.dihedrals_data = dihedrals_data 

self.imdihedrals_data = imdihedrals_data 

 

def __str__(self): 

""" 

returns a string of lammps data input file 

""" 

lines = [] 

# title 

lines.append("Data file generated by rubicon\n") 

 

# count 

lines.append("{} atoms".format(self.natoms)) 

lines.append("{} bonds".format(self.nbonds)) 

lines.append("{} angles".format(self.nangles)) 

if self.ndih > 0: 

lines.append("{} dihedrals".format(self.ndih)) 

if self.nimdihs > 0: 

lines.append("{} impropers".format(self.nimdihs)) 

 

# types 

lines.append("\n{} atom types".format(self.natom_types)) 

lines.append("{} bond types".format(self.nbond_types)) 

lines.append("{} angle types".format(self.nangle_types)) 

if self.ndih > 0: 

lines.append("{} dihedral types".format(self.ndih_types)) 

if self.nimdihs > 0: 

lines.append("{} improper types".format(self.nimdih_types)) 

 

# box size 

lines.append("\n{} {} xlo xhi\n{} {} ylo yhi\n{} {} zlo zhi".format( 

self.box_size[0][0], self.box_size[0][1], 

self.box_size[1][0], self.box_size[1][1], 

self.box_size[2][0], self.box_size[2][1])) 

 

# masses 

self.set_lines_from_list(lines, "Masses", self.atomic_masses) 

 

# coefficients 

self.set_lines_from_list(lines, "Pair Coeffs", self.pair_coeffs) 

self.set_lines_from_list(lines, "Bond Coeffs", self.bond_coeffs) 

self.set_lines_from_list(lines, "Angle Coeffs", self.angle_coeffs) 

if self.ndih > 0: 

self.set_lines_from_list(lines, "Dihedral Coeffs", 

self.dihedral_coeffs) 

if self.nimdihs > 0: 

self.set_lines_from_list(lines, "Improper Coeffs", 

self.improper_coeffs) 

 

# data 

self.set_lines_from_list(lines, "Atoms", self.atoms_data) 

self.set_lines_from_list(lines, "Bonds", self.bonds_data) 

self.set_lines_from_list(lines, "Angles", self.angles_data) 

if self.ndih > 0: 

self.set_lines_from_list(lines, "Dihedrals", self.dihedrals_data) 

if self.nimdihs > 0: 

self.set_lines_from_list(lines, "Impropers", self.imdihedrals_data) 

return '\n'.join(lines) 

 

@staticmethod 

def get_param_coeff(forcefield, param_name): 

""" 

get the parameter coefficients and mapping from the force field. 

 

Args: 

forcefield (ForceField): ForceField object 

param_name (string): name of the parameter for which 

the coefficients are to be set. 

 

Returns: 

[[parameter id, value1, value2, ... ], ... ] and 

{parameter key: parameter id, ...} 

""" 

if hasattr(forcefield, param_name): 

param = getattr(forcefield, param_name) 

if param: 

param_coeffs = [[i + 1] + list(v) for i, v in 

enumerate(param.values())] 

param_map = dict([(k, i) for i, k in enumerate(param.keys())]) 

return param_coeffs, param_map 

else: 

return [], {} 

else: 

raise AttributeError 

 

@staticmethod 

def get_atoms_data(mols, mols_number, molecule, atomic_masses_dict, 

topologies): 

""" 

Return the atoms data. 

 

Args: 

mols (list): list of Molecule objects. 

mols_number (list): number of each type of molecule in mols list. 

molecule (Molecule): the molecule assembled from the molecules 

in the mols list. 

topologies (list): list of Topology objects, one for each molecule 

type in mols list 

 

Returns: 

atoms_data: [[atom id, mol type, atom type, charge, x, y, z], ... ] 

molid_to_atomid: [ [global atom id 1, id 2, ..], ...], the 

index will be the global mol id 

""" 

atom_to_mol = {} 

molid_to_atomid = [] 

atoms_data = [] 

nmols = len(mols) 

shift = [0] + [len(mols[i]) * mols_number[i] for i in range(nmols - 1)] 

# set up map atom_id --> [mol_type, local atom id in the mol] in mols 

# set up map gobal molecule id --> [[atom_id,...],...] 

for mol_type in range(nmols): 

natoms = len(mols[mol_type]) 

for num_mol_id in range(mols_number[mol_type]): 

tmp = [] 

for mol_atom_id in range(natoms): 

atom_id = num_mol_id * natoms + mol_atom_id + shift[ 

mol_type] 

atom_to_mol[atom_id] = [mol_type, mol_atom_id] 

tmp.append(atom_id) 

molid_to_atomid.append(tmp) 

# set atoms data from the final molecule assembly(the packed molecule 

# obtained from packmol using the molecules from mols list and the 

# molecule numbers from mol_number list). 

# atom id, mol id, atom type, charge from topology, x, y, z 

for i, site in enumerate(molecule): 

atom_type = atomic_masses_dict[site.specie.symbol][0] 

# atom_type = molecule.symbol_set.index(site.species_string) + 1 

atom_id = i + 1 

mol_type = atom_to_mol[i][0] + 1 

mol_atom_id = atom_to_mol[i][1] + 1 

charge = 0.0 

if hasattr(topologies[0], "charges"): 

if topologies[mol_type - 1].charges: 

charge = topologies[mol_type - 1].charges[mol_atom_id - 1][ 

0] 

atoms_data.append([atom_id, mol_type, atom_type, charge, 

site.x, site.y, site.z]) 

return atoms_data, molid_to_atomid 

 

@staticmethod 

def get_param_data(param_name, param_map, mols, mols_number, topologies, 

molid_to_atomid): 

""" 

set the data for the parameter named param_name from the topology. 

 

Args: 

param_name (string): parameter name, example: "bonds" 

param_map (dict): 

{ mol_type: {parameter_key : unique parameter id, ... }, ... } 

example: {0: {("c1","c2"): 1}} ==> c1-c2 bond in mol_type=0 

has the global id of 1 

mols (list): list of molecules. 

mols_number (list): number of each type of molecule in mols list. 

topologies (list): list of Topology objects, one for each molecule 

type in mols list 

molid_to_atomid (list): [ [gloabal atom id 1, id 2, ..], ...], 

the index is the global mol id 

 

Returns: 

[ [parameter id, parameter type, global atom id1, global atom id2, ...], ... ] 

""" 

if hasattr(topologies[0], param_name) and getattr(topologies[0], 

param_name): 

nmols = len(mols) 

shift = [0] + [len(mols[i]) * mols_number[i] for i in 

range(nmols - 1)] 

mol_id = 0 

# set the map param_to_mol: 

# {global param_id :[global mol id, mol_type, local param id in the param], ... } 

param_to_mol = {} 

for mol_type in range(nmols): 

param_obj = getattr(topologies[mol_type], param_name) 

nparams = len(param_obj) 

for num_mol_id in range(mols_number[mol_type]): 

mol_id += 1 

for mol_param_id in range(nparams): 

param_id = num_mol_id * nparams + mol_param_id + shift[ 

mol_type] 

param_to_mol[param_id] = [mol_id - 1, mol_type, 

mol_param_id] 

# set the parameter data using the topology info 

# example: loop over all bonds in the system 

param_data = [] 

for param_id, pinfo in param_to_mol.items(): 

mol_id = pinfo[0] # global molecule id 

mol_type = pinfo[1] # type of molecule 

mol_param_id = pinfo[2] # local parameter id in that molecule 

# example: get the bonds list for mol_type molecule 

param_obj = getattr(topologies[mol_type], param_name) 

# connectivity info(local atom ids and type) for the parameter with the local id 

# 'mol_param_id'. example: single bond = [i, j, bond_type] 

param = param_obj[mol_param_id] 

param_atomids = [] 

# loop over local atom ids that constitute the parameter 

# for the molecule type, mol_type 

# example: single bond = [i,j,bond_label] 

for atomid in param[:-1]: 

# local atom id to global atom id 

global_atom_id = molid_to_atomid[mol_id][atomid - 1] 

param_atomids.append(global_atom_id + 1) 

param_type = param[-1] 

# example: get the unique number id for the bond_type 

param_type_id = param_map[param_type] 

param_data.append( 

[param_id + 1, param_type_id + 1] + param_atomids) 

return param_data 

else: 

return [] 

 

@staticmethod 

def from_forcefield_and_topology(mols, mols_number, box_size, molecule, 

forcefield, 

topologies): 

""" 

Return LammpsForceFieldData object from force field and topology info. 

 

Args: 

mols (list): List of Molecule objects 

mols_number (list): List of number of molecules of each 

molecule type in mols 

box_size (list): [[x_min,x_max], [y_min,y_max], [z_min,z_max]] 

molecule (Molecule): The molecule that is assembled from mols 

and mols_number 

forcefield (ForceFiled): Force filed information 

topologies (list): List of Topology objects, one for each 

molecule type in mols. 

 

Returns: 

LammpsForceFieldData 

""" 

# set the coefficients and map from the force field 

bond_coeffs, bond_map = LammpsForceFieldData.get_param_coeff( 

forcefield, "bonds") 

angle_coeffs, angle_map = LammpsForceFieldData.get_param_coeff( 

forcefield, "angles") 

pair_coeffs, _ = LammpsForceFieldData.get_param_coeff(forcefield, 

"vdws") 

dihedral_coeffs, dihedral_map = LammpsForceFieldData.get_param_coeff( 

forcefield, "dihedrals") 

improper_coeffs, imdihedral_map = LammpsForceFieldData.get_param_coeff( 

forcefield, "imdihedrals") 

# atoms data, topology used for setting charge if present 

molecule = get_boxed_molecule(molecule, box_size) 

natoms, natom_types, box_size, atomic_masses_dict = LammpsData.get_basic_system_info( 

molecule) 

atoms_data, molid_to_atomid = LammpsForceFieldData.get_atoms_data(mols, 

mols_number, 

molecule, 

atomic_masses_dict, 

topologies) 

# set the other data from the molecular topologies 

bonds_data = LammpsForceFieldData.get_param_data("bonds", bond_map, 

mols, mols_number, 

topologies, 

molid_to_atomid) 

angles_data = LammpsForceFieldData.get_param_data("angles", angle_map, 

mols, mols_number, 

topologies, 

molid_to_atomid) 

dihedrals_data = LammpsForceFieldData.get_param_data("dihedrals", 

dihedral_map, 

mols, 

mols_number, 

topologies, 

molid_to_atomid) 

imdihedrals_data = LammpsForceFieldData.get_param_data("imdihedrals", 

imdihedral_map, 

mols, 

mols_number, 

topologies, 

molid_to_atomid) 

return LammpsForceFieldData(box_size, atomic_masses_dict.values(), 

pair_coeffs, bond_coeffs, 

angle_coeffs, dihedral_coeffs, 

improper_coeffs, atoms_data, 

bonds_data, angles_data, dihedrals_data, 

imdihedrals_data) 

 

@staticmethod 

def from_file(data_file): 

""" 

Return LammpsForceFieldData object from the data file 

 

Args: 

data_file (string): the data file name 

 

Returns: 

LammpsForceFieldData 

""" 

atomic_masses = [] # atom_type(starts from 1): mass 

box_size = [] 

pair_coeffs = [] 

bond_coeffs = [] 

angle_coeffs = [] 

dihedral_coeffs = [] 

improper_coeffs = [] 

atoms_data = [] 

bonds_data = [] 

angles_data = [] 

dihedral_data = [] 

imdihedral_data = [] 

types_pattern = re.compile("^\s*(\d+)\s+([a-zA-Z]+)\s+types$") 

# atom_id, mol_id, atom_type, charge, x, y, z 

atoms_pattern = re.compile( 

"^\s*(\d+)\s+(\d+)\s+(\d+)\s+([0-9eE\.+-]+)\s+(" 

"[0-9eE\.+-]+)\s+([0-9eE\.+-]+)\s+([0-9eE\.+-]+)\w*") 

masses_pattern = re.compile("^\s*(\d+)\s+([0-9\.]+)$") 

box_pattern = re.compile( 

"^\s*([0-9eE\.+-]+)\s+([0-9eE\.+-]+)\s+[xyz]lo\s+[xyz]hi") 

# id, value1, value2 

general_coeff_pattern = re.compile( 

"^\s*(\d+)\s+([0-9\.]+)\s+([0-9\.]+)$") 

# id, type, atom_id1, atom_id2 

bond_data_pattern = re.compile("^\s*(\d+)\s+(\d+)\s+(\d+)\s+(\d+)\s*$") 

# id, type, atom_id1, atom_id2, atom_id3 

angle_data_pattern = re.compile( 

"^\s*(\d+)\s+(\d+)\s+(\d+)\s+(\d+)\s+(\d+)\s*$") 

# id, type, atom_id1, atom_id2, atom_id3, atom_id4 

dihedral_data_pattern = re.compile( 

"^\s*(\d+)\s+(\d+)\s+(\d+)\s+(\d+)\s+(\d+)\s+(\d+)\s*$") 

read_pair_coeffs = False 

read_bond_coeffs = False 

read_angle_coeffs = False 

read_dihedral_coeffs = False 

read_improper_coeffs = False 

with open(data_file) as df: 

for line in df: 

if types_pattern.search(line): 

m = types_pattern.search(line) 

if m.group(2) == "atom": 

natom_types = int(m.group(1)) 

if m.group(2) == "bond": 

nbond_types = int(m.group(1)) 

if m.group(2) == "angle": 

nangle_types = int(m.group(1)) 

if m.group(2) == "dihedral": 

ndihedral_types = int(m.group(1)) 

if m.group(2) == "improper": 

nimproper_types = int(m.group(1)) 

if masses_pattern.search(line): 

m = masses_pattern.search(line) 

atomic_masses.append([int(m.group(1)), float(m.group(2))]) 

if box_pattern.search(line): 

m = box_pattern.search(line) 

box_size.append([float(m.group(1)), float(m.group(2))]) 

if "Pair Coeffs" in line: 

read_pair_coeffs = True 

continue 

if read_pair_coeffs: 

m = general_coeff_pattern.search(line) 

if m: 

pair_coeffs.append([int(m.group(1)), float(m.group(2)), 

float(m.group(3))]) 

read_pair_coeffs = False if len( 

pair_coeffs) >= natom_types else True 

if "Bond Coeffs" in line: 

read_bond_coeffs = True 

continue 

if read_bond_coeffs: 

m = general_coeff_pattern.search(line) 

if m: 

bond_coeffs.append([int(m.group(1)), float(m.group(2)), 

float(m.group(3))]) 

read_bond_coeffs = False if len( 

bond_coeffs) >= nbond_types else True 

if "Angle Coeffs" in line: 

read_angle_coeffs = True 

continue 

if read_angle_coeffs: 

m = general_coeff_pattern.search(line) 

if m: 

angle_coeffs.append( 

[int(m.group(1)), float(m.group(2)), 

float(m.group(3))]) 

read_angle_coeffs = False if len( 

angle_coeffs) >= nangle_types else True 

if "Dihedral Coeffs" in line: 

read_dihedral_coeffs = True 

continue 

if read_dihedral_coeffs: 

m = general_coeff_pattern.search(line) 

if m: 

dihedral_coeffs.append( 

[int(m.group(1))] + [float(i) for i in 

m.groups()[1:]]) 

read_dihedral_coeffs = False if len( 

dihedral_coeffs) >= ndihedral_types \ 

else True 

if "Improper Coeffs" in line: 

read_improper_coeffs = True 

continue 

if read_improper_coeffs: 

m = general_coeff_pattern.search(line) 

if m: 

improper_coeffs.append( 

[int(m.group(1))] + [float(i) for i in 

m.groups()[1:]]) 

read_improper_coeffs = False if len( 

improper_coeffs) > nimproper_types else True 

if atoms_pattern.search(line): 

m = atoms_pattern.search(line) 

# atom id, mol id, atom type 

line_data = [int(i) for i in m.groups()[:3]] 

# charge, x, y, z, vx, vy, vz ... 

line_data.extend([float(i) for i in m.groups()[3:]]) 

atoms_data.append(line_data) 

if bond_data_pattern.search(line): 

m = bond_data_pattern.search(line) 

bonds_data.append([int(i) for i in m.groups()]) 

if angle_data_pattern.search(line): 

m = angle_data_pattern.search(line) 

angles_data.append([int(i) for i in m.groups()]) 

if dihedral_data_pattern.search(line): 

m = dihedral_data_pattern.search(line) 

dihedral_data.append([int(i) for i in m.groups()]) 

return LammpsForceFieldData(box_size, atomic_masses, pair_coeffs, 

bond_coeffs, angle_coeffs, 

dihedral_coeffs, improper_coeffs, 

atoms_data, bonds_data, angles_data, 

dihedral_data, imdihedral_data) 

 

 

def get_boxed_molecule(input_structure, box_size): 

""" 

Return a boxed molecule. 

 

Args: 

input_structure(Molecule/Structure): either Molecule of Structure object 

box_size (list): [[x_min, x_max], [y_min, y_max], [z_min, z_max]] 

 

Returns: 

Structure 

""" 

box_lengths = [min_max[1] - min_max[0] for min_max in box_size] 

# be defensive about the box size 

if isinstance(input_structure, Molecule): 

boxed_molecule = input_structure.get_boxed_structure(*box_lengths) 

elif isinstance(input_structure, Structure): 

max_length = max(input_structure.lattice.abc) 

max_box_size = max(box_lengths) 

boxed_molecule = Molecule.from_sites(input_structure.sites) 

if max_length < max_box_size: 

boxed_molecule = \ 

boxed_molecule.get_boxed_structure(*box_lengths) 

else: 

boxed_molecule = \ 

boxed_molecule.get_boxed_structure(max_length, 

max_length, max_length) 

else: 

raise ValueError("molecule must be an object of Molecule or " 

"Structure ") 

return boxed_molecule