/
molecule.py
1914 lines (1642 loc) · 77.4 KB
/
molecule.py
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
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#
# @BEGIN LICENSE
#
# Psi4: an open-source quantum chemistry software package
#
# Copyright (c) 2007-2019 The Psi4 Developers.
#
# The copyrights for code used from other parties are included in
# the corresponding files.
#
# This file is part of Psi4.
#
# Psi4 is free software; you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as published by
# the Free Software Foundation, version 3.
#
# Psi4 is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License along
# with Psi4; if not, write to the Free Software Foundation, Inc.,
# 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
#
# @END LICENSE
#
import os
import hashlib
import collections
import numpy as np
import qcelemental as qcel
from .util import parse_dertype
from .libmintsmolecule import *
from .testing import compare_values, compare_integers, compare_molrecs
from .bfs import BFS
class Molecule(LibmintsMolecule):
"""Class to store the elements, coordinates, fragmentation pattern,
charge, multiplicity of a molecule. Largely replicates psi4's libmints
Molecule class, developed by Justin M. Turney and Andy M. Simmonett
with incremental improvements by other psi4 developers. Major
This class extends `qcdb.LibmintsMolecule` and occasionally
`psi4.core.Molecule` itself.
"""
def __init__(self,
molinit=None,
dtype=None,
geom=None,
elea=None,
elez=None,
elem=None,
mass=None,
real=None,
elbl=None,
name=None,
units='Angstrom',
input_units_to_au=None,
fix_com=None,
fix_orientation=None,
fix_symmetry=None,
fragment_separators=None,
fragment_charges=None,
fragment_multiplicities=None,
molecular_charge=None,
molecular_multiplicity=None,
comment=None,
provenance=None,
connectivity=None,
enable_qm=True,
enable_efp=True,
missing_enabled_return_qm='none',
missing_enabled_return_efp='none',
missing_enabled_return='error',
tooclose=0.1,
zero_ghost_fragments=False,
nonphysical=False,
mtol=1.e-3,
verbose=1):
"""Initialize Molecule object from LibmintsMolecule"""
super(Molecule, self).__init__()
if molinit is not None or geom is not None:
if isinstance(molinit, dict):
molrec = molinit
elif isinstance(molinit, str):
compound_molrec = qcel.molparse.from_string(
molstr=molinit,
dtype=dtype,
name=name,
fix_com=fix_com,
fix_orientation=fix_orientation,
fix_symmetry=fix_symmetry,
return_processed=False,
enable_qm=enable_qm,
enable_efp=enable_efp,
missing_enabled_return_qm=missing_enabled_return_qm,
missing_enabled_return_efp=missing_enabled_return_efp,
verbose=verbose)
molrec = compound_molrec['qm']
elif molinit is None and geom is not None:
molrec = qcel.molparse.from_arrays(
geom=geom,
elea=elea,
elez=elez,
elem=elem,
mass=mass,
real=real,
elbl=elbl,
name=name,
units=units,
input_units_to_au=input_units_to_au,
fix_com=fix_com,
fix_orientation=fix_orientation,
fix_symmetry=fix_symmetry,
fragment_separators=fragment_separators,
fragment_charges=fragment_charges,
fragment_multiplicities=fragment_multiplicities,
molecular_charge=molecular_charge,
molecular_multiplicity=molecular_multiplicity,
comment=comment,
provenance=provenance,
connectivity=connectivity,
domain='qm',
missing_enabled_return=missing_enabled_return,
tooclose=tooclose,
zero_ghost_fragments=zero_ghost_fragments,
nonphysical=nonphysical,
mtol=mtol,
verbose=verbose)
# ok, got the molrec dictionary; now build the thing
self._internal_from_dict(molrec, verbose=verbose)
# The comment line
self.tagline = ""
def __str__(self):
text = """ ==> qcdb Molecule %s <==\n\n""" % (self.name())
text += """ => %s <=\n\n""" % (self.tagline)
text += self.create_psi4_string_from_molecule()
return text
def __setattr__(self, name, value):
"""Function to overload setting attributes to allow geometry
variable assigment as if member data.
"""
if 'all_variables' in self.__dict__:
if name.upper() in self.__dict__['all_variables']:
self.set_variable(name, value)
super(Molecule, self).__setattr__(name, value)
def __getattr__(self, name):
"""Function to overload accessing attribute contents to allow
retrival of geometry variable values as if member data.
"""
if 'all_variables' in self.__dict__ and name.upper() in self.__dict__['all_variables']:
return self.get_variable(name)
else:
raise AttributeError
@classmethod
def init_with_xyz(cls, xyzfilename, no_com=False, no_reorient=False, contentsNotFilename=False):
"""Pull information from an XYZ file. No fragment info detected.
Bohr/Angstrom pulled from first line if available. Charge,
multiplicity, tagline pulled from second line if available. Body
accepts atom symbol or atom charge in first column. Arguments
*no_com* and *no_reorient* can be used to turn off shift and
rotation. If *xyzfilename* is a string of the contents of an XYZ
file, rather than the name of a file, set *contentsNotFilename*
to ``True``.
>>> H2O = qcdb.Molecule.init_with_xyz('h2o.xyz')
"""
raise FeatureDeprecated(
"""qcdb.Molecule.init_with_xyz. Replace with: qcdb.Molecule.from_string(..., dtype='xyz+')""")
@classmethod
def init_with_mol2(cls, xyzfilename, no_com=False, no_reorient=False, contentsNotFilename=False):
"""Pull information from a MOl2 file. No fragment info detected.
Bohr/Angstrom pulled from first line if available. Charge,
multiplicity, tagline pulled from second line if available. Body
accepts atom symbol or atom charge in first column. Arguments
*no_com* and *no_reorient* can be used to turn off shift and
rotation. If *xyzfilename* is a string of the contents of an XYZ
file, rather than the name of a file, set *contentsNotFilename*
to ``True``.
NOTE: chg/mult NYI
>>> H2O = qcdb.Molecule.init_with_mol2('h2o.mol2')
"""
instance = cls()
instance.lock_frame = False
instance.PYmove_to_com = not no_com
instance.PYfix_orientation = no_reorient
if contentsNotFilename:
text = xyzfilename.splitlines()
else:
try:
infile = open(xyzfilename, 'r')
except IOError:
raise ValidationError(
"""Molecule::init_with_mol2: given filename '%s' does not exist.""" % (xyzfilename))
if os.stat(xyzfilename).st_size == 0:
raise ValidationError("""Molecule::init_with_mol2: given filename '%s' is blank.""" % (xyzfilename))
text = infile.readlines()
# fixed-width regex ((?=[ ]*-?\d+)[ -\d]{5})
v2000 = re.compile(r'^((?=[ ]*\d+)[ \d]{3})((?=[ ]*\d+)[ \d]{3})(.*)V2000\s*$')
vend = re.compile(r'^\s*M\s+END\s*$')
NUMBER = "((?:[-+]?\\d*\\.\\d+(?:[DdEe][-+]?\\d+)?)|(?:[-+]?\\d+\\.\\d*(?:[DdEe][-+]?\\d+)?))"
xyzM = re.compile(
r'^(?:\s*)' + NUMBER + r'(?:\s+)' + NUMBER + r'(?:\s+)' + NUMBER + r'(?:\s+)([A-Z](?:[a-z])?)(?:\s+)(.*)',
re.IGNORECASE)
## now charge and multiplicity
# $chargem = 0 ; $multm = 1 ;
#while (<MOL>) {
#if (/CHARGE/) { $chargem = <MOL> ; chop($chargem) ;}
#if (/MULTIPLICITY/) { $multm = <MOL> ; chop($multm) }
# } # end while charge and multiplicity
if not text:
raise ValidationError("Molecule::init_with_mol2: file blank")
# Try to match header/footer
if vend.match(text[-1]):
pass
else:
raise ValidationError("Molecule::init_with_mol2: Malformed file termination\n%s" % (text[-1]))
sysname = '_'.join(text[0].strip().split())
comment = text[2].strip()
if comment:
instance.tagline = sysname + ' ' + comment
else:
instance.tagline = sysname
#instance.tagline = text[0].strip() + ' ' + text[2].strip()
fileUnits = 'Angstrom' # defined for MOL
#instance.set_molecular_charge(int(xyz2.match(text[1]).group(1)))
#instance.set_multiplicity(int(xyz2.match(text[1]).group(2)))
if v2000.match(text[3]):
fileNatom = int(v2000.match(text[3]).group(1))
fileNbond = int(v2000.match(text[3]).group(2))
else:
raise ValidationError("Molecule::init_with_mol2: Malformed fourth line\n%s" % (text[3]))
if fileNatom < 1:
raise ValidationError("Molecule::init_with_mol2: Malformed Natom\n%s" % (str(fileNatom)))
# Next line begins the useful information.
for i in range(fileNatom):
try:
if xyzM.match(text[4 + i]):
fileX = float(xyzM.match(text[4 + i]).group(1))
fileY = float(xyzM.match(text[4 + i]).group(2))
fileZ = float(xyzM.match(text[4 + i]).group(3))
fileAtom = xyzM.match(text[4 + i]).group(4).upper()
# Check that the atom symbol is valid
z = qcel.periodictable.to_Z(fileAtom)
# Add it to the molecule.
instance.add_atom(z, fileX, fileY, fileZ, fileAtom, qcel.periodictable.to_mass(fileAtom), z)
else:
raise ValidationError("Molecule::init_with_mol2: Malformed atom information line %d." % (i + 5))
except IndexError:
raise ValidationError(
"Molecule::init_with_mol2: Expected atom in file at line %d.\n%s" % (i + 5, text[i + 4]))
# We need to make 1 fragment with all atoms
instance.fragments.append([0, fileNatom - 1])
instance.fragment_types.append('Real')
instance.fragment_charges.append(instance.molecular_charge())
instance.fragment_multiplicities.append(instance.multiplicity())
# Set the units properly
instance.PYunits = fileUnits
if fileUnits == 'Bohr':
instance.PYinput_units_to_au = 1.0
elif fileUnits == 'Angstrom':
instance.PYinput_units_to_au = 1.0 / qcel.constants.bohr2angstroms
instance.update_geometry()
return instance
def save_string_xyz(self, save_ghosts=True, save_natom=False):
"""Save a string for a XYZ-style file.
>>> H2OH2O.save_string_xyz()
6
-2 3 water_dimer
O -1.551007000000 -0.114520000000 0.000000000000
H -1.934259000000 0.762503000000 0.000000000000
H -0.599677000000 0.040712000000 0.000000000000
O 1.350625000000 0.111469000000 0.000000000000
H 1.680398000000 -0.373741000000 -0.758561000000
H 1.680398000000 -0.373741000000 0.758561000000
"""
factor = 1.0 if self.PYunits == 'Angstrom' else qcel.constants.bohr2angstroms
N = self.natom()
if not save_ghosts:
N = 0
for i in range(self.natom()):
if self.Z(i):
N += 1
text = ''
if save_natom:
text += "%d\n" % (N)
text += '%d %d %s\n' % (self.molecular_charge(), self.multiplicity(), self.tagline)
for i in range(self.natom()):
[x, y, z] = self.atoms[i].compute()
if save_ghosts or self.Z(i):
text += '%2s %17.12f %17.12f %17.12f\n' % ((self.symbol(i) if self.Z(i) else "Gh"), \
x * factor, y * factor, z * factor)
return text
def save_xyz(self, filename, save_ghosts=True, save_natom=True):
"""Save an XYZ file.
>>> H2OH2O.save_xyz('h2o.xyz')
"""
outfile = open(filename, 'w')
outfile.write(self.save_string_xyz(save_ghosts, save_natom))
outfile.close()
def format_molecule_for_numpy(self, npobj=True):
"""Returns a NumPy array of the non-dummy atoms of the geometry
in Cartesian coordinates in Angstroms with element encoded as
atomic number. If *npobj* is False, returns representation of
NumPy array.
"""
factor = 1.0 if self.PYunits == 'Angstrom' else qcel.constants.bohr2angstroms
self.update_geometry()
# TODO fn title is format_mol... but return args not compatible
geo = []
for i in range(self.natom()):
[x, y, z] = self.atoms[i].compute()
geo.append([self.Z(i), x * factor, y * factor, z * factor])
nparr = np.array(geo)
return nparr if npobj else np.array_repr(nparr)
def format_molecule_for_psi4(self):
"""Returns string of molecule definition block."""
text = 'molecule mol {\n'
for line in self.create_psi4_string_from_molecule().splitlines():
text += ' ' + line + '\n'
text += '}\n'
return text
def format_molecule_for_qchem_old(self, mixedbas=True):
"""Returns geometry section of input file formatted for Q-Chem.
For ghost atoms, prints **Gh** as elemental symbol, with expectation
that element identity will be established in mixed basis section.
For ghost atoms when *mixedbas* is False, prints @ plus element symbol.
prints whole dimer for unCP mono when called dir (as opposed to passing thru str
no frag markers
"""
factor = 1.0 if self.PYunits == 'Angstrom' else qcel.constants.bohr2angstroms
text = ""
text += '$molecule\n'
text += '%d %d\n' % (self.molecular_charge(), self.multiplicity())
for i in range(self.natom()):
[x, y, z] = self.atoms[i].compute()
if mixedbas:
text += '%2s ' % (self.symbol(i) if self.Z(i) else "Gh")
else:
text += '%-3s ' % (('' if self.Z(i) else '@') + self.symbol(i))
text += '%17.12f %17.12f %17.12f\n' % (x * factor, y * factor, z * factor)
text += '$end\n\n'
# prepare molecule keywords to be set as c-side keywords
options = collections.defaultdict(lambda: collections.defaultdict(dict))
#options['QCHEM'['QCHEM_CHARGE']['value'] = self.molecular_charge()
#options['QCHEM'['QCHEM_MULTIPLICITY']['value'] = self.multiplicity()
options['QCHEM']['QCHEM_INPUT_BOHR']['value'] = False
#options['QCHEM']['QCHEM_COORDINATES']['value'] = 'CARTESIAN'
#SYM_IGNORE equiv to no_reorient, no_com, symmetry c1
options['QCHEM']['QCHEM_INPUT_BOHR']['clobber'] = True
return text, options
def format_molecule_for_psi4_xyz(self):
"""not much examined
"""
text = ""
if self.nallatom():
factor = 1.0 if self.PYunits == 'Angstrom' else qcel.constants.bohr2angstroms
# append units and any other non-default molecule keywords
text += "units Angstrom\n"
#text += " units %-s\n" % ("Angstrom" if self.units() == 'Angstrom' else "Bohr")
if not self.PYmove_to_com:
text += "no_com\n"
if self.PYfix_orientation:
text += "no_reorient\n"
# append atoms and coordentries and fragment separators with charge and multiplicity
Pfr = 0
for fr in range(self.nfragments()):
if self.fragment_types[fr] == 'Absent' and not self.has_zmatrix():
continue
text += "%s%s%d %d\n" % ("" if Pfr == 0 else "--\n", "#" if self.fragment_types[fr] == 'Ghost'
or self.fragment_types[fr] == 'Absent' else "", self.fragment_charges[fr],
self.fragment_multiplicities[fr])
Pfr += 1
for at in range(self.fragments[fr][0], self.fragments[fr][1] + 1):
if self.fragment_types[fr] == 'Absent' or self.fsymbol(at) == "X":
pass
else:
if self.fZ(at):
text += "%-8s" % (self.flabel(at))
else:
text += "%-8s" % ("Gh(" + self.flabel(at) + ")")
[x, y, z] = self.full_atoms[at].compute()
text += '%17.12f %17.12f %17.12f\n' % \
(x * factor, y * factor, z * factor)
text += "\n"
wtext = 'molecule mol {\n'
for line in text.splitlines():
wtext += ' ' + line + '\n'
wtext += '}\n'
return wtext
def format_molecule_for_molpro(self):
"""
"""
factor = 1.0 if self.PYunits == 'Angstrom' else qcel.constants.bohr2angstroms
# TODO keep fix_or? # Jan 2015 turning off fix_or
#self.fix_orientation(True)
#self.PYmove_to_com = False
self.update_geometry()
text = ""
text += 'angstrom\n'
text += 'geometry={\n'
dummy = []
for i in range(self.natom()):
[x, y, z] = self.atoms[i].compute()
text += '%-2s %17.12f %17.12f %17.12f\n' % (self.symbol(i), \
x * factor, y * factor, z * factor)
if not self.Z(i):
dummy.append(str(i + 1)) # Molpro atom number is 1-indexed
text += '}\n\n'
text += 'SET,CHARGE=%d\n' % (self.molecular_charge())
text += 'SET,SPIN=%d\n' % (self.multiplicity() - 1) # Molpro wants (mult-1)
if len(dummy) > 0:
text += 'dummy,' + ','.join(dummy) + '\n'
return text
def format_molecule_for_cfour(self):
"""Function to print Molecule in a form readable by Cfour.
"""
self.update_geometry()
factor = 1.0 if self.PYunits == 'Angstrom' else qcel.constants.bohr2angstroms
#factor = 1.0 if self.PYunits == 'Bohr' else 1.0/psi_bohr2angstroms
text = 'auto-generated by qcdb from molecule %s\n' % (self.tagline)
# append atoms and coordentries
for i in range(self.natom()):
[x, y, z] = self.atoms[i].compute()
text += '%-2s %17.12f %17.12f %17.12f\n' % ((self.symbol(i) if self.Z(i) else "GH"), \
x * factor, y * factor, z * factor)
#for fr in range(self.nfragments()):
# if self.fragment_types[fr] == 'Absent':
# pass
# else:
# for at in range(self.fragments[fr][0], self.fragments[fr][1] + 1):
# [x, y, z] = self.atoms[at].compute()
# text += '%-2s %17.12f %17.12f %17.12f\n' % ((self.symbol(at) if self.Z(at) else "GH"), \
# x * factor, y * factor, z * factor)
text += '\n'
# prepare molecule keywords to be set as c-side keywords
options = collections.defaultdict(lambda: collections.defaultdict(dict))
options['CFOUR']['CFOUR_CHARGE']['value'] = self.molecular_charge()
options['CFOUR']['CFOUR_MULTIPLICITY']['value'] = self.multiplicity()
options['CFOUR']['CFOUR_UNITS']['value'] = 'ANGSTROM'
#options['CFOUR']['CFOUR_UNITS']['value'] = 'BOHR'
options['CFOUR']['CFOUR_COORDINATES']['value'] = 'CARTESIAN'
#options['CFOUR']['CFOUR_SUBGROUP']['value'] = self.symmetry_from_input().upper()
#print self.inertia_tensor()
#print self.inertial_system()
options['CFOUR']['CFOUR_CHARGE']['clobber'] = True
options['CFOUR']['CFOUR_MULTIPLICITY']['clobber'] = True
options['CFOUR']['CFOUR_UNITS']['clobber'] = True
options['CFOUR']['CFOUR_COORDINATES']['clobber'] = True
return text, options
def format_basis_for_cfour(self, puream):
"""Function to print the BASIS=SPECIAL block for Cfour according
to the active atoms in Molecule. Special short basis names
are used by Psi4 libmints GENBAS-writer in accordance with
Cfour constraints.
"""
text = ''
cr = 1
for fr in range(self.nfragments()):
if self.fragment_types[fr] == 'Absent':
pass
else:
for at in range(self.fragments[fr][0], self.fragments[fr][1] + 1):
text += """%s:P4_%d\n""" % (self.symbol(at).upper(), cr)
cr += 1
text += '\n'
options = collections.defaultdict(lambda: collections.defaultdict(dict))
options['CFOUR']['CFOUR_BASIS']['value'] = 'SPECIAL'
options['CFOUR']['CFOUR_SPHERICAL']['value'] = puream
options['CFOUR']['CFOUR_BASIS']['clobber'] = True
options['CFOUR']['CFOUR_SPHERICAL']['clobber'] = True
options['CFOUR']['CFOUR_BASIS']['superclobber'] = True
options['CFOUR']['CFOUR_SPHERICAL']['superclobber'] = True
return text, options
def format_molecule_for_orca(self):
"""
Format the molecule into an orca xyz format
"""
options = collections.defaultdict(lambda: collections.defaultdict(dict))
self.update_geometry()
factor = 1.0 if self.PYunits == 'Angstrom' else qcel.constants.bohr2angstroms
text = ""
text += '* xyz {} {}\n'.format(self.molecular_charge(), self.multiplicity())
n_frags = self.nfragments()
for fr in range(n_frags):
if self.fragment_types[fr] == 'Absent':
pass
else:
for at in range(self.fragments[fr][0], self.fragments[fr][1] + 1):
if self.fragment_types[fr] == 'Ghost':
# TODO: add support for ghost atoms
# atom += ':'
continue
x, y, z = self.atoms[at].compute()
atom = self.symbol(at)
if n_frags > 1:
text += ' {:2s}({:d}) {:> 17.12f} {:> 17.12f} {:> 17.12f}\n'.format(\
atom, fr + 1, x * factor, y * factor, z * factor)
else:
text += ' {:2s} {:> 17.12f} {:> 17.12f} {:> 17.12f}\n'.format(\
atom, x * factor, y * factor, z * factor)
text += '*'
return text, options
def format_molecule_for_qchem(self, mixedbas=True):
"""Returns geometry section of input file formatted for Q-Chem.
For ghost atoms, prints **Gh** as elemental symbol, with expectation
that element identity will be established in mixed basis section.
For ghost atoms when *mixedbas* is False, prints @ plus element symbol.
candidate modeled after psi4_xyz so that absent fragments observed force xyz
"""
text = ""
if self.nallatom():
factor = 1.0 if self.PYunits == 'Angstrom' else qcel.constants.bohr2angstroms
Pfr = 0
# any general starting notation here <<<
text += '$molecule\n'
text += '%d %d\n' % (self.molecular_charge(), self.multiplicity())
# >>>
for fr in range(self.nfragments()):
if self.fragment_types[fr] == 'Absent' and not self.has_zmatrix():
continue
# any fragment marker here <<<
if self.nactive_fragments() > 1:
# this only distinguishes Real frags so Real/Ghost don't get
# fragmentation. may need to change
text += """--\n"""
# >>>
# any fragment chgmult here <<<
if self.nactive_fragments() > 1:
text += """{}{} {}\n""".format('!' if self.fragment_types[fr] in ['Ghost', 'Absent'] else '',
self.fragment_charges[fr], self.fragment_multiplicities[fr])
# >>>
Pfr += 1
for at in range(self.fragments[fr][0], self.fragments[fr][1] + 1):
if self.fragment_types[fr] == 'Absent' or self.fsymbol(at) == "X":
pass
else:
if self.fZ(at):
# label for real live atom <<<
text += """{:>3s} """.format(self.fsymbol(at))
# >>>
else:
# label for ghost atom <<<
text += """{:>3s} """.format('Gh' if mixedbas else ('@' + self.fsymbol(at)))
# >>>
[x, y, z] = self.full_atoms[at].compute()
# Cartesian coordinates <<<
text += """{:>17.12f} {:>17.12f} {:>17.12f}\n""".format(x * factor, y * factor, z * factor)
# >>>
# any general finishing notation here <<<
text += '$end\n\n'
# >>>
# prepare molecule keywords to be set as c-side keywords
options = collections.defaultdict(lambda: collections.defaultdict(dict))
#options['QCHEM'['QCHEM_CHARGE']['value'] = self.molecular_charge()
#options['QCHEM'['QCHEM_MULTIPLICITY']['value'] = self.multiplicity()
options['QCHEM']['QCHEM_INPUT_BOHR']['value'] = False
#options['QCHEM']['QCHEM_COORDINATES']['value'] = 'CARTESIAN'
if (not self.PYmove_to_com) or self.PYfix_orientation:
options['QCHEM']['QCHEM_SYM_IGNORE']['value'] = True
#SYM_IGNORE equiv to no_reorient, no_com, symmetry c1
options['QCHEM']['QCHEM_INPUT_BOHR']['clobber'] = True
options['QCHEM']['QCHEM_SYM_IGNORE']['clobber'] = True
return text, options
def format_molecule_for_cfour_old(self):
"""Function to print Molecule in a form readable by Cfour. This
version works as long as zmat is composed entirely of variables,
not internal values, while cartesian is all internal values,
no variables. Cutting off this line of development because,
with getting molecules after passing through libmints Molecule,
all zmats with dummies (Cfour's favorite kind) have already been
converted into cartesian. Next step, if this line was pursued
would be to shift any zmat internal values to external and any
cartesian external values to internal.
"""
text = ''
text += 'auto-generated by qcdb from molecule %s\n' % (self.tagline)
## append units and any other non-default molecule keywords
#text += " units %-s\n" % ("Angstrom" if self.units() == 'Angstrom' else "Bohr")
#if not self.PYmove_to_com:
# text += " no_com\n"
#if self.PYfix_orientation:
# text += " no_reorient\n"
# append atoms and coordentries and fragment separators with charge and multiplicity
Pfr = 0
isZMat = False
isCart = False
for fr in range(self.nfragments()):
if self.fragment_types[fr] == 'Absent' and not self.has_zmatrix():
continue
# text += "%s %s%d %d\n" % (
# "" if Pfr == 0 else " --\n",
# "#" if self.fragment_types[fr] == 'Ghost' or self.fragment_types[fr] == 'Absent' else "",
# self.fragment_charges[fr], self.fragment_multiplicities[fr])
Pfr += 1
for at in range(self.fragments[fr][0], self.fragments[fr][1] + 1):
if type(self.full_atoms[at]) == ZMatrixEntry:
isZMat = True
elif type(self.full_atoms[at]) == CartesianEntry:
isCart = True
if self.fragment_types[fr] == 'Absent':
text += "%s" % ("X")
elif self.fZ(at) or self.fsymbol(at) == "X":
text += "%s" % (self.fsymbol(at))
else:
text += "%s" % ("GH") # atom info is lost + self.fsymbol(at) + ")")
text += "%s" % (self.full_atoms[at].print_in_input_format_cfour())
text += "\n"
# append any coordinate variables
if len(self.geometry_variables):
for vb, val in self.geometry_variables.items():
text += """%s=%.10f\n""" % (vb, val)
text += "\n"
# prepare molecule keywords to be set as c-side keywords
options = collections.defaultdict(lambda: collections.defaultdict(dict))
options['CFOUR']['CFOUR_CHARGE']['value'] = self.molecular_charge()
options['CFOUR']['CFOUR_MULTIPLICITY']['value'] = self.multiplicity()
options['CFOUR']['CFOUR_UNITS']['value'] = self.units()
if isZMat and not isCart:
options['CFOUR']['CFOUR_COORDINATES']['value'] = 'INTERNAL'
elif isCart and not isZMat:
options['CFOUR']['CFOUR_COORDINATES']['value'] = 'CARTESIAN'
else:
raise ValidationError("""Strange mix of Cartesian and ZMatrixEntries in molecule unsuitable for Cfour.""")
return text, options
def format_molecule_for_nwchem(self):
"""
"""
factor = 1.0 if self.PYunits == 'Angstrom' else qcel.constants.bohr2angstroms
text = ""
text += '%d %d %s\n' % (self.molecular_charge(), self.multiplicity(), self.tagline)
for i in range(self.natom()):
[x, y, z] = self.atoms[i].compute()
text += '%4s %17.12f %17.12f %17.12f\n' % (("" if self.Z(i) else 'Bq') + self.symbol(i), \
x * factor, y * factor, z * factor)
return text
pass
# if symm print M2OUT "nosym\nnoorient\n";
# print DIOUT "angstrom\ngeometry={\n";
def inertia_tensor(self, masswt=True, zero=ZERO):
"""Compute inertia tensor.
>>> print H2OH2O.inertia_tensor()
[[8.704574864178731, -8.828375721817082, 0.0], [-8.828375721817082, 280.82861714077666, 0.0], [0.0, 0.0, 281.249500988553]]
"""
return self.inertia_tensor_partial(range(self.natom()), masswt, zero)
def inertia_tensor_partial(self, part, masswt=True, zero=ZERO):
"""Compute inertia tensor based on atoms in *part*.
"""
tensor = [[0, 0, 0], [0, 0, 0], [0, 0, 0]]
for i in part:
if masswt:
# I(alpha, alpha)
tensor[0][0] += self.mass(i) * (self.y(i) * self.y(i) + self.z(i) * self.z(i))
tensor[1][1] += self.mass(i) * (self.x(i) * self.x(i) + self.z(i) * self.z(i))
tensor[2][2] += self.mass(i) * (self.x(i) * self.x(i) + self.y(i) * self.y(i))
# I(alpha, beta)
tensor[0][1] -= self.mass(i) * self.x(i) * self.y(i)
tensor[0][2] -= self.mass(i) * self.x(i) * self.z(i)
tensor[1][2] -= self.mass(i) * self.y(i) * self.z(i)
else:
# I(alpha, alpha)
tensor[0][0] += self.y(i) * self.y(i) + self.z(i) * self.z(i)
tensor[1][1] += self.x(i) * self.x(i) + self.z(i) * self.z(i)
tensor[2][2] += self.x(i) * self.x(i) + self.y(i) * self.y(i)
# I(alpha, beta)
tensor[0][1] -= self.x(i) * self.y(i)
tensor[0][2] -= self.x(i) * self.z(i)
tensor[1][2] -= self.y(i) * self.z(i)
# mirror
tensor[1][0] = tensor[0][1]
tensor[2][0] = tensor[0][2]
tensor[2][1] = tensor[1][2]
# Check the elements for zero and make them a hard zero.
for i in range(3):
for j in range(3):
if math.fabs(tensor[i][j]) < zero:
tensor[i][j] = 0.0
return tensor
def inertial_system_partial(self, part, masswt=True, zero=ZERO):
"""Solve inertial system based on atoms in *part*"""
return diagonalize3x3symmat(self.inertia_tensor_partial(part, masswt, zero))
def inertial_system(self, masswt=True, zero=ZERO):
"""Solve inertial system"""
return diagonalize3x3symmat(self.inertia_tensor(masswt, zero))
def print_ring_planes(self, entity1, entity2, entity3=None, entity4=None):
"""(reals only, 1-indexed)
"""
pass
# TODO allow handle lines
text = ""
summ = []
#for entity in [entity1, entity2, entity3, entity4]:
for item in [entity1, entity2]:
text += """\n ==> Entity %s <==\n\n""" % (item)
# convert plain atoms into list and move from 1-indexed to 0-indexed
entity = []
try:
for idx in item:
entity.append(idx - 1)
except TypeError:
entity = [item - 1]
if len(entity) == 1:
dim = 'point'
elif len(entity) == 2:
dim = 'line'
else:
dim = 'plane'
# compute centroid
cent = [0.0, 0.0, 0.0]
for at in entity:
cent = add(cent, self.xyz(at))
cent = scale(cent, 1.0 / len(entity))
text += ' Centroid: %14.8f %14.8f %14.8f [Angstrom]\n' % \
(cent[0] * qcel.constants.bohr2angstroms, \
cent[1] * qcel.constants.bohr2angstroms, \
cent[2] * qcel.constants.bohr2angstroms)
text += ' Centroid: %14.8f %14.8f %14.8f [Bohr]\n' % \
(cent[0], cent[1], cent[2])
if dim == 'point':
summ.append({'dim': dim, 'geo': cent, 'cent': cent})
# TODO: figure out if should be using mass-weighted
self.translate(scale(cent, -1))
evals, evecs = self.inertial_system_partial(entity, masswt=False)
midx = evals.index(max(evals))
text += ' Normal Vector: %14.8f %14.8f %14.8f [unit]\n' % \
(evecs[0][midx], evecs[1][midx], evecs[2][midx])
text += ' Normal Vector: %14.8f %14.8f %14.8f [unit]\n' % \
(evecs[0][midx] + cent[0], evecs[1][midx] + cent[1], evecs[2][midx] + cent[2])
xplane = [evecs[0][midx], evecs[1][midx], evecs[2][midx], \
-1.0 * (evecs[0][midx] * cent[0] + evecs[1][midx] * cent[1] + evecs[2][midx] * cent[2])]
text += ' Eqn. of Plane: %14.8f %14.8f %14.8f %14.8f [Ai + Bj + Ck + D = 0]\n' % \
(xplane[0], xplane[1], xplane[2], xplane[3])
dtemp = math.sqrt(evecs[0][midx] * evecs[0][midx] + evecs[1][midx] * evecs[1][midx] +
evecs[2][midx] * evecs[2][midx])
hessplane = [evecs[0][midx] / dtemp, evecs[1][midx] / dtemp, evecs[2][midx] / dtemp, xplane[3] / dtemp]
hessplane2 = [xplane[0] / dtemp, xplane[1] / dtemp, xplane[2] / dtemp, xplane[3] / dtemp]
text += ' Eqn. of Plane: %14.8f %14.8f %14.8f %14.8f [Ai + Bj + Ck + D = 0] H\n' % \
(hessplane[0], hessplane[1], hessplane[2], hessplane[3])
text += ' Eqn. of Plane: %14.8f %14.8f %14.8f %14.8f [Ai + Bj + Ck + D = 0] H2\n' % \
(hessplane2[0], hessplane2[1], hessplane2[2], hessplane2[3])
self.translate(cent)
if dim == 'plane':
summ.append({'dim': dim, 'geo': xplane, 'cent': cent})
#print summ
text += """\n ==> 1 (%s) vs. 2 (%s) <==\n\n""" % (summ[0]['dim'], summ[1]['dim'])
#if summ[0]['dim'] == 'plane' and summ[1]['dim'] == 'point':
# cent = summ[1]['geo']
# plane = summ[0]['geo']
# print cent, plane
#
# D = math.fabs(plane[0] * cent[0] + plane[1] * cent[1] + plane[2] * cent[2] + plane[3]) / \
# math.sqrt(plane[0] * plane[0] + plane[1] * plane[1] + plane[2] * plane[2])
# text += ' Pt to Plane: %14.8f [Angstrom]\n' % (D * psi_bohr2angstroms)
#if summ[0]['dim'] == 'plane' and summ[1]['dim'] == 'plane':
if summ[0]['dim'] == 'plane' and (summ[1]['dim'] == 'plane' or summ[1]['dim'] == 'point'):
cent1 = summ[0]['cent']
cent2 = summ[1]['cent']
plane1 = summ[0]['geo']
#plane2 = summ[1]['geo']
distCC = distance(cent1, cent2)
text += ' Distance from Center of %s to Center of %s: %14.8f [Angstrom]\n' % \
('2', '1', distCC * qcel.constants.bohr2angstroms)
distCP = math.fabs(plane1[0] * cent2[0] + plane1[1] * cent2[1] + plane1[2] * cent2[2] + plane1[3])
# distCP expression has a denominator that's one since plane constructed from unit vector
text += ' Distance from Center of %s to Plane of %s: %14.8f [Angstrom]\n' % \
('2', '1', distCP * qcel.constants.bohr2angstroms)
distCPC = math.sqrt(distCC * distCC - distCP * distCP)
text += ' Distance from Center of %s to Center of %s along Plane of %s: %14.8f [Angstrom]\n' % \
('2', '1', '1', distCPC * qcel.constants.bohr2angstroms)
print(text)
# text = " Interatomic Distances (Angstroms)\n\n"
# for i in range(self.natom()):
# for j in range(i + 1, self.natom()):
# eij = sub(self.xyz(j), self.xyz(i))
# dist = norm(eij) * psi_bohr2angstroms
# text += " Distance %d to %d %-8.3lf\n" % (i + 1, j + 1, dist)
# text += "\n\n"
# return text
def rotor_type(self, tol=FULL_PG_TOL):
"""Returns the rotor type.
>>> H2OH2O.rotor_type()
RT_ASYMMETRIC_TOP
"""
evals, evecs = diagonalize3x3symmat(self.inertia_tensor())
evals = sorted(evals)
rot_const = [
1.0 / evals[0] if evals[0] > 1.0e-6 else 0.0, 1.0 / evals[1] if evals[1] > 1.0e-6 else 0.0,
1.0 / evals[2] if evals[2] > 1.0e-6 else 0.0
]
# Determine degeneracy of rotational constants.
degen = 0
for i in range(2):
for j in range(i + 1, 3):
if degen >= 2:
continue
rabs = math.fabs(rot_const[i] - rot_const[j])
tmp = rot_const[i] if rot_const[i] > rot_const[j] else rot_const[j]
if rabs > ZERO:
rel = rabs / tmp
else:
rel = 0.0
if rel < tol:
degen += 1
#print "\tDegeneracy is %d\n" % (degen)
# Determine rotor type
if self.natom() == 1:
rotor_type = 'RT_ATOM'
elif rot_const[0] == 0.0:
rotor_type = 'RT_LINEAR' # 0 < IB == IC inf > B == C
elif degen == 2:
rotor_type = 'RT_SPHERICAL_TOP' # IA == IB == IC A == B == C
elif degen == 1:
if (rot_const[1] - rot_const[2]) < 1.0e-6:
rotor_type = 'RT_PROLATE_SYMMETRIC_TOP' # IA < IB == IC A > B == C
elif (rot_const[0] - rot_const[1]) < 1.0e-6:
rotor_type = 'RT_OBLATE_SYMMETRIC_TOP' # IA == IB < IC A == B > C
else:
rotor_type = 'RT_ASYMMETRIC_TOP' # IA < IB < IC A > B > C
return rotor_type
def center_of_charge(self):
"""Computes center of charge of molecule (does not translate molecule).
>>> H2OH2O.center_of_charge()
[-0.073339893272065401, 0.002959783555632145, 0.0]
"""
ret = [0.0, 0.0, 0.0]
total_c = 0.0
for at in range(self.natom()):
c = self.charge(at)
ret = add(ret, scale(self.xyz(at), c))
total_c += c
ret = scale(ret, 1.0 / total_c)
return ret
def move_to_coc(self):
"""Moves molecule to center of charge
"""
coc = scale(self.center_of_charge(), -1.0)
self.translate(coc)
def rotational_symmetry_number(self):
"""Number of unique orientations of the rigid molecule that only interchange identical atoms.
Notes
-----
Source http://cccbdb.nist.gov/thermo.asp (search "symmetry number")
"""
pg = self.get_full_point_group()
pg = self.full_point_group_with_n()
if pg in ['ATOM', 'C1', 'Ci', 'Cs', 'C_inf_v']:
sigma = 1
elif pg == 'D_inf_h':