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
|
# -*- coding: utf-8 -*-
#
# Copyright (c) 2023, the cclib development team
#
# This file is part of cclib (http://cclib.github.io) and is distributed under
# the terms of the BSD 3-Clause License.
"""Test single point logfiles in cclib."""
import datetime
import os
import unittest
import numpy
import packaging
from common import get_minimum_carbon_separation
from skip import skipForParser
from skip import skipForLogfile
__filedir__ = os.path.realpath(os.path.dirname(__file__))
class GenericSPTest(unittest.TestCase):
"""Generic restricted single point unittest"""
# Molecular mass of DVB in mD, and expected precision.
molecularmass = 130078.25
mass_precision = 0.10
# In STO-3G, H has 1, C has 5 (1 S and 4 SP).
nbasisdict = {1:1, 6:5}
# Approximate B3LYP energy of dvb after SCF in STO-3G.
b3lyp_energy = -10365
# Approximate energy of the innermost molecular orbital of DVB with
# B3LYP/STO-3G (from Q-Chem 5.4).
b3lyp_moenergy = -272.60365543
b3lyp_moenergy_delta = 7.55e-2
# Overlap first two atomic orbitals.
overlap01 = 0.24
# Generally, one criteria for SCF energy convergence.
num_scf_criteria = 1
def testnatom(self):
"""Is the number of atoms equal to 20?"""
assert self.data.natom == 20
def testatomnos(self):
"""Are the atomnos correct?"""
# The nuclear charges should be integer values in a NumPy array.
assert numpy.alltrue([numpy.issubdtype(atomno, numpy.signedinteger)
for atomno in self.data.atomnos])
assert self.data.atomnos.dtype.char == 'i'
assert self.data.atomnos.shape == (20,)
assert sum(self.data.atomnos == 6) + sum(self.data.atomnos == 1) == 20
@skipForParser('DALTON', 'DALTON has a very low accuracy for the printed values of all populations (2 decimals rounded in a weird way), so let it slide for now')
@skipForParser('FChk', 'The parser is still being developed so we skip this test')
@skipForParser('GAMESSDAT', 'We are not going to implement atom charges in this version.')
@skipForLogfile('Jaguar/basicJaguar7', 'We did not print the atomic partial charges in the unit tests for this version')
@skipForLogfile('Molpro/basicMolpro2006', "These tests were run a long time ago and since we don't have access to Molpro 2006 anymore, we can skip this test (it is tested in 2012)")
@skipForParser('Turbomole','The parser is still being developed so we skip this test')
def testatomcharges(self):
"""Are atomic charges consistent with natom?"""
for atomcharge_type in self.data.atomcharges:
charges = self.data.atomcharges[atomcharge_type]
natom = self.data.natom
assert len(charges) == natom, f"len(atomcharges['{atomcharge_type}']) = {len(charges)}, natom = {natom}"
@skipForParser('DALTON', 'DALTON has a very low accuracy for the printed values of all populations (2 decimals rounded in a weird way), so let it slide for now')
@skipForLogfile("FChk/basicQChem5.2", "not printed for Q-Chem")
@skipForLogfile("FChk/basicQChem5.4", "not printed for Q-Chem")
@skipForParser('GAMESSDAT', 'We are not sure about the specific type of atom charges, it is best to skip the test for now.')
@skipForLogfile('Jaguar/basicJaguar7', 'We did not print the atomic partial charges in the unit tests for this version')
@skipForLogfile('Molpro/basicMolpro2006', "These tests were run a long time ago and since we don't have access to Molpro 2006 anymore, we can skip this test (it is tested in 2012)")
@skipForParser('Turbomole','The parser is still being developed so we skip this test')
def testatomcharges_mulliken(self):
"""Do Mulliken atomic charges sum to zero?"""
charges = self.data.atomcharges["mulliken"]
assert abs(sum(charges)) < 1.0e-2
@skipForParser('ADF', 'Lowdin charges not present by default')
@skipForParser('DALTON', 'DALTON has a very low accuracy for the printed values of all populations (2 decimals rounded in a weird way), so let it slide for now')
@skipForParser('FChk', 'The parser is still being developed so we skip this test')
@skipForParser('GAMESSDAT', 'We are not sure about the specific type of atom charges, it is best to skip the test for now.')
@skipForParser('Gaussian', 'Lowdin charges not present by default')
@skipForParser('Jaguar', 'Lowdin charges not present by default')
@skipForParser('NWChem', 'Lowdin charges not present by default')
@skipForParser('Molcas', 'Lowdin charges not present by default')
@skipForParser('Molpro', 'Lowdin charges not present by default')
@skipForParser('QChem', 'Lowdin charges not present by default')
@skipForParser('Turbomole','The parser is still being developed so we skip this test')
def testatomcharges_lowdin(self):
"""Do Lowdin atomic charges sum to zero?"""
charges = self.data.atomcharges["lowdin"]
assert abs(sum(charges)) < 1.0e-2
@skipForParser('ADF', 'Hirshfeld charges not implemented')
@skipForParser('DALTON', 'DALTON has a very low accuracy for the printed values of all populations (2 decimals rounded in a weird way), so let it slide for now')
@skipForParser('FChk', 'The parser is still being developed so we skip this test')
@skipForParser('Gaussian', 'Hirshfeld charges not implemented')
@skipForParser('GAMESS', 'Hirshfeld charges not implemented')
@skipForParser('GAMESSUK', 'Hirshfeld charges not implemented')
@skipForParser('GAMESSDAT', 'We are not sure about the specific type of atom charges, it is best to skip the test for now.')
@skipForParser('Jaguar', 'Hirshfeld charges not implemented')
@skipForParser('NWChem', 'Hirshfeld charges not implemented')
@skipForParser('Molcas', 'Hirshfeld charges not implemented')
@skipForParser('Molpro', 'Hirshfeld charges not implemented')
@skipForLogfile('ORCA/basicORCA4.1', 'This needs to be moved to regressions')
@skipForParser('Psi4', 'Hirshfeld charges not implemented')
@skipForParser('QChem', 'Hirshfeld charges not implemented')
@skipForParser('Turbomole','The parser is still being developed so we skip this test')
def testatomcharges_hirshfeld(self):
"""Do Hirshfeld atomic charges sum to roughly zero?"""
charges = self.data.atomcharges["hirshfeld"]
assert abs(sum(charges)) < 4.0e-3
def testatomcoords(self):
"""Are the dimensions of atomcoords 1 x natom x 3?"""
expected_shape = (1, self.data.natom, 3)
assert self.data.atomcoords.shape == expected_shape
@skipForParser('GAMESSDAT', 'Vectors need some calculations to transform them. Current mm value is 2.54')
def testatomcoords_units(self):
"""Are atomcoords consistent with Angstroms?"""
min_carbon_dist = get_minimum_carbon_separation(self.data)
dev = abs(min_carbon_dist - 1.34)
assert dev < 0.03, f"Minimum carbon dist is {min_carbon_dist:.2f} (not 1.34)"
@skipForParser('GAMESSDAT', 'Neither charge nor mult exists in the files.')
@skipForParser('Molcas', 'missing mult')
def testcharge_and_mult(self):
"""Are the charge and multiplicity correct?"""
assert self.data.charge == 0
assert self.data.mult == 1
def testnbasis(self):
"""Is the number of basis set functions correct?"""
count = sum([self.nbasisdict[n] for n in self.data.atomnos])
assert self.data.nbasis == count
@skipForParser('ADF', 'ADF parser does not extract atombasis')
@skipForLogfile('Jaguar/basicJaguar7', 'Data file does not contain enough information. Can we make a new one?')
@skipForParser('Molcas','The parser is still being developed so we skip this test')
@skipForParser('Turbomole','The parser is still being developed so we skip this test')
def testatombasis(self):
"""Are the indices in atombasis the right amount and unique?"""
all = []
for i, atom in enumerate(self.data.atombasis):
assert len(atom) == self.nbasisdict[self.data.atomnos[i]]
all += atom
# Test if there are as many indices as atomic orbitals.
assert len(all) == self.data.nbasis
# Check if all are different (every orbital indexed once).
assert len(set(all)) == len(all)
@skipForLogfile("FChk/basicQChem5.2", "Q-Chem doesn't print SCF energy to fchk")
@skipForLogfile("FChk/basicQChem5.4", "Q-Chem doesn't print SCF energy to fchk")
@skipForParser('GAMESS', 'atommasses not implemented yet')
@skipForParser('GAMESSUK', 'atommasses not implemented yet')
@skipForParser('GAMESSDAT', 'Atommasses implemented, but it does not pass the test.')
@skipForParser('Jaguar', 'atommasses not implemented yet')
@skipForParser('Molcas','The parser is still being developed so we skip this test')
@skipForParser('Molpro', 'atommasses not implemented yet')
@skipForLogfile('Psi4/basicPsi4.0b5', 'atommasses not implemented yet')
@skipForParser('QChem', 'atommasses not implemented yet')
@skipForParser('Turbomole','The parser is still being developed so we skip this test')
def testatommasses(self):
"""Do the atom masses sum up to the molecular mass?"""
mm = 1000*sum(self.data.atommasses)
msg = f"Molecule mass: {mm:f} not {self.molecularmass:f} +- {self.mass_precision:f}mD"
assert abs(mm-self.molecularmass) < self.mass_precision, msg
@skipForParser('Turbomole','The parser is still being developed so we skip this test')
def testcoreelectrons(self):
"""Are the coreelectrons all 0?"""
ans = numpy.zeros(self.data.natom, 'i')
numpy.testing.assert_array_equal(self.data.coreelectrons, ans)
@skipForParser('FChk', 'Formatted checkpoint files do not have a section for symmetry')
@skipForParser('GAMESSDAT', 'Mosyms do not exist in the file')
@skipForParser('Molcas','The parser is still being developed so we skip this test')
@skipForParser('Molpro', '?')
def testsymlabels(self):
"""Are all the symmetry labels either Ag/u or Bg/u?"""
sumwronglabels = sum([x not in ['Ag', 'Bu', 'Au', 'Bg'] for x in self.data.mosyms[0]])
assert sumwronglabels == 0
def testhomos(self):
"""Is the index of the HOMO equal to 34?"""
numpy.testing.assert_array_equal(
self.data.homos,
numpy.array([34], "i"),
f"{numpy.array_repr(self.data.homos)} != array([34],'i')",
)
@skipForParser('FChk', 'Formatted Checkpoint files do not have a section for SCF energy')
@skipForParser('GAMESSDAT', 'Scfvalues probably do not exist in the file')
def testscfvaluetype(self):
"""Are scfvalues and its elements the right type??"""
assert isinstance(self.data.scfvalues, list)
assert isinstance(self.data.scfvalues[0], numpy.ndarray)
@skipForLogfile("FChk/basicQChem5.2", "Q-Chem doesn't print SCF energy to fchk")
@skipForLogfile("FChk/basicQChem5.4", "Q-Chem doesn't print SCF energy to fchk")
@skipForParser('GAMESSDAT', 'Scfenergies probably do not exist in the file')
def testscfenergy(self):
"""Is the SCF energy within the target?"""
assert abs(self.data.scfenergies[-1]-self.b3lyp_energy) < 40
@skipForParser('FChk', 'Formatted Checkpoint files do not have a section for SCF convergence')
@skipForParser('GAMESSDAT', 'Scftargets probably do not exist in the file')
def testscftargetdim(self):
"""Do the scf targets have the right dimensions?"""
assert self.data.scftargets.shape == (len(self.data.scfvalues), len(self.data.scfvalues[0][0]))
@skipForParser('FChk', 'Formatted Checkpoint files do not have a section for SCF convergence')
@skipForParser('GAMESSDAT', 'Scftargets probably do not exist in the file')
def testscftargets(self):
"""Are correct number of SCF convergence criteria being parsed?"""
assert len(self.data.scftargets[0]) == self.num_scf_criteria
def testlengthmoenergies(self):
"""Is the number of evalues equal to nmo?"""
if hasattr(self.data, "moenergies"):
assert len(self.data.moenergies[0]) == self.data.nmo
def testtypemoenergies(self):
"""Is moenergies a list containing one numpy array?"""
if hasattr(self.data, "moenergies"):
assert isinstance(self.data.moenergies, list)
assert isinstance(self.data.moenergies[0], numpy.ndarray)
@skipForParser('GAMESSDAT', 'Moenergies probably do not exist in the file')
@skipForLogfile('Gaussian/basicGaussian16/dvb_sp_no.out', 'no energies for natural orbitals')
@skipForLogfile('Turbomole/basicTurbomole5.9/dvb_sp_symm', 'delta of 7.4, everything else ok')
def testfirstmoenergy(self):
"""Is the lowest energy molecular orbital within the target?"""
assert abs(self.data.moenergies[0][0]-self.b3lyp_moenergy) < self.b3lyp_moenergy_delta
@skipForParser('DALTON', 'mocoeffs not implemented yet')
@skipForLogfile('Jaguar/basicJaguar7', 'Data file does not contain enough information. Can we make a new one?')
@skipForParser('Turbomole', 'Use of symmetry has reduced the number of mo coeffs')
def testdimmocoeffs(self):
"""Are the dimensions of mocoeffs equal to 1 x nmo x nbasis?"""
if hasattr(self.data, "mocoeffs"):
assert isinstance(self.data.mocoeffs, list)
assert len(self.data.mocoeffs) == 1
assert self.data.mocoeffs[0].shape == (self.data.nmo, self.data.nbasis)
@skipForParser('DALTON', 'mocoeffs not implemented yet')
@skipForLogfile('Jaguar/basicJaguar7', 'Data file does not contain enough information. Can we make a new one?')
def testfornoormo(self):
"""Do we have NOs or MOs?"""
assert hasattr(self.data, "nocoeffs") or hasattr(self.data, "mocoeffs")
def testdimnoccnos(self):
"""Is the length of nooccnos equal to nmo?"""
if hasattr(self.data, "nooccnos"):
assert isinstance(self.data.nooccnos, numpy.ndarray)
assert len(self.data.nooccnos) == self.data.nmo
def testdimnocoeffs(self):
"""Are the dimensions of nocoeffs equal to nmo x nmo?"""
if hasattr(self.data, "nocoeffs"):
assert isinstance(self.data.nocoeffs, numpy.ndarray)
assert self.data.nocoeffs.shape == (self.data.nmo, self.data.nmo)
@skipForParser('DALTON', 'To print: **INTEGRALS\n.PROPRI')
@skipForLogfile('FChk/basicGaussian09', 'Only available in QChem')
@skipForLogfile('FChk/basicGaussian16', 'Only available in QChem')
@skipForParser('GAMESSDAT','Aooverlaps probably do not exist in the file.')
@skipForParser('Molcas','The parser is still being developed so we skip this test')
@skipForParser('Psi4', 'Psi4 does not currently have the option to print the overlap matrix')
@skipForParser('QChem', 'QChem cannot print the overlap matrix')
@skipForParser('Turbomole','The parser is still being developed so we skip this test')
def testaooverlaps(self):
"""Are the dims and values of the overlap matrix correct?"""
assert self.data.aooverlaps.shape == (self.data.nbasis, self.data.nbasis)
# The matrix is symmetric.
row = self.data.aooverlaps[0,:]
col = self.data.aooverlaps[:,0]
assert sum(col - row) == 0.0
# All values on diagonal should be exactly one.
for i in range(self.data.nbasis):
assert self.data.aooverlaps[i,i] == 1.0
# Check some additional values that don't seem to move around between programs.
assert abs(self.data.aooverlaps[0, 1] - self.overlap01) < 0.01
assert abs(self.data.aooverlaps[1, 0] - self.overlap01) < 0.01
assert round(abs(self.data.aooverlaps[3, 0]), 7) == 0
assert round(abs(self.data.aooverlaps[0, 3]), 7) == 0
def testoptdone(self):
"""There should be no optdone attribute set."""
assert not hasattr(self.data, 'optdone')
@skipForParser('ADF', 'Not implemented yes')
@skipForParser('DALTON', 'Not implemented yes')
@skipForParser('FChk', 'Rotational constants are never written to fchk files')
@skipForParser('GAMESS', 'Not implemented yes')
@skipForParser('GAMESSUK', 'Not implemented yet')
@skipForParser('GAMESSDAT', 'Not implemented yet')
@skipForParser('Jaguar', 'Not implemented yet')
@skipForParser('Molcas', 'Not implemented yes')
@skipForParser('Molpro', 'Not implemented yes')
@skipForParser('NWChem', 'Not implemented yes')
@skipForParser('ORCA', 'Not implemented yes')
@skipForParser('Psi4', 'Not implemented yes')
@skipForParser('QChem', 'Not implemented yes')
@skipForParser('Turbomole', 'Not implemented yes')
def testrotconsts(self):
"""A single geometry leads to single set of rotational constants."""
assert self.data.rotconsts.shape == (1, 3)
# taken from Gaussian16/dvb_sp.out
ref = [4.6266363, 0.6849065, 0.5965900]
numpy.testing.assert_allclose(self.data.rotconsts[0], ref, rtol=0, atol=1.0e-3)
@skipForParser('FChk', 'The parser is still being developed so we skip this test')
@skipForParser('Gaussian', 'Logfile needs to be updated')
@skipForParser('Jaguar', 'No dipole moments in the logfile')
@skipForParser('Molcas','The parser is still being developed so we skip this test')
def testmoments(self):
"""Does the dipole and possible higher molecular moments look reasonable?"""
# The reference point is always a vector, but not necessarily the
# origin or center of mass. In this case, however, the center of mass
# is at the origin, so we now what to expect.
reference = self.data.moments[0]
assert len(reference) == 3
for x in reference:
assert x == 0.0
# Length and value of dipole moment should always be correct (zero for this test).
dipole = self.data.moments[1]
assert len(dipole) == 3
for d in dipole:
assert round(abs(d), 7) == 0
# If the quadrupole is there, we can expect roughly -50B for the XX moment,
# -50B for the YY moment and and -60B for the ZZ moment.
if len(self.data.moments) > 2:
quadrupole = self.data.moments[2]
assert len(quadrupole) == 6
assert abs(quadrupole[0] - -50) < 2.5
assert abs(quadrupole[3] - -50) < 2.5
assert abs(quadrupole[5] - -60) < 3
# If the octupole is there, it should have 10 components and be zero.
if len(self.data.moments) > 3:
octupole = self.data.moments[3]
assert len(octupole) == 10
for m in octupole:
assert abs(m) < 0.001
# The hexadecapole should have 15 elements, an XXXX component of around -1900 Debye*ang^2,
# a YYYY component of -330B and a ZZZZ component of -50B.
if len(self.data.moments) > 4:
hexadecapole = self.data.moments[4]
assert len(hexadecapole) == 15
assert abs(hexadecapole[0] - -1900) < 90
assert abs(hexadecapole[10] - -330) < 11
assert abs(hexadecapole[14] - -50) < 2.5
# The are 21 unique 32-pole moments, and all are zero in this test case.
if len(self.data.moments) > 5:
moment32 = self.data.moments[5]
assert len(moment32) == 21
for m in moment32:
assert m == 0.0
@skipForParser('ADF', 'reading basis set names is not implemented')
@skipForParser('GAMESSDAT', 'Basis set not implemented in this version.')
@skipForParser('GAMESSUK', 'reading basis set names is not implemented')
@skipForParser('Molcas', 'reading basis set names is not implemented')
@skipForParser('Psi4', 'reading basis set names is not implemented')
def testmetadata_basis_set(self):
"""Does metadata have expected keys and values?"""
assert self.data.metadata["basis_set"].lower() == "sto-3g"
@skipForParser('ADF', 'reading input file contents and name is not implemented')
@skipForParser('DALTON', 'reading input file contents and name is not implemented')
@skipForParser('FChk', 'Formatted checkpoint files do not have an input file section')
@skipForParser('GAMESS', 'reading input file contents and name is not implemented')
@skipForParser('GAMESSUK', 'reading input file contents and name is not implemented')
@skipForParser('GAMESSDAT', 'reading input file contents and name is not implemented')
@skipForParser('Gaussian', 'reading input file contents and name is not implemented')
@skipForParser('Jaguar', 'reading input file contents and name is not implemented')
@skipForParser('Molcas', 'reading input file contents and name is not implemented')
@skipForParser('Molpro', 'reading input file contents and name is not implemented')
@skipForParser('NWChem', 'reading input file contents and name is not implemented')
@skipForParser('Psi4', 'reading input file contents and name is not implemented')
@skipForParser('QChem', 'reading input file contents and name is not implemented')
@skipForParser('Turbomole', 'reading input file contents and name is not implemented')
def testmetadata_input_file(self):
"""Does metadata have expected keys and values?"""
assert "input_file_contents" in self.data.metadata
# TODO make input file names consistent where possible, though some
# programs do not allow arbitrary file extensions; for example, DALTON
# must end in `dal`.
assert "dvb_sp.in" in self.data.metadata["input_file_name"]
def testmetadata_methods(self):
"""Does metadata have expected keys and values?"""
# TODO implement and unify across parsers; current values are [],
# ["HF"], ["RHF"], and ["DFT"]
assert "methods" in self.data.metadata
def testmetadata_package(self):
"""Does metadata have expected keys and values?"""
# TODO How can the value be tested when the package name comes from
# the parser and isn't stored on ccData?
assert "package" in self.data.metadata
@skipForParser('FChk', 'Formatted Checkpoint files do not have section for legacy package version')
@skipForParser('GAMESSDAT', 'Files do not contain information about the legacy package version')
def testmetadata_legacy_package_version(self):
"""Does metadata have expected keys and values?"""
# TODO Test specific values for each unit test.
assert "legacy_package_version" in self.data.metadata
@skipForParser('FChk', 'Formatted Checkpoint files do not have section for package version')
@skipForParser('GAMESSDAT', 'Files do not contain information about the package version')
def testmetadata_package_version(self):
"""Does metadata have expected keys and values?"""
# TODO Test specific values for each unit test.
assert isinstance(packaging.version.parse(self.data.metadata["package_version"]), packaging.version.Version)
@skipForParser('FChk', 'point group symmetry cannot be printed')
@skipForParser('GAMESSDAT', 'Files probably do not contain information about symmetry_detected')
@skipForParser('Molcas', 'reading point group symmetry and name is not implemented')
@skipForParser('Molpro', 'reading point group symmetry and name is not implemented')
@skipForParser('Turbomole', 'reading point group symmetry and name is not implemented')
def testmetadata_symmetry_detected(self):
"""Does metadata have expected keys and values?"""
assert self.data.metadata["symmetry_detected"] == "c2h"
@skipForParser('FChk', 'point group symmetry cannot be printed')
@skipForParser('GAMESSDAT', 'Files probably do not contain information about symmetry_used')
@skipForParser('Molcas', 'reading point group symmetry and name is not implemented')
@skipForParser('Molpro', 'reading point group symmetry and name is not implemented')
@skipForParser('Turbomole', 'reading point group symmetry and name is not implemented')
def testmetadata_symmetry_used(self):
"""Does metadata have expected keys and values?"""
assert self.data.metadata["symmetry_used"] == "c2h"
@skipForParser('ADF', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('DALTON', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('FChk', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('GAMESS', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('GAMESSUK', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('GAMESSUS', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('GAMESSDAT', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('Jaguar', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('Molcas', ' reading cpu/wall time is not implemented for this parser')
@skipForParser('Molpro', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('NWChem', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('Psi3', 'reading cpu/wall time is not implemented for this parser')
@skipForParser('Psi4', 'reading cpu/wall time is not implemented for this parser')
def testmetadata_times(self):
"""Does metadata have expected keys and values of correct types?"""
if "wall_time" in self.data.metadata:
assert self.data.metadata["wall_time"]
assert all(isinstance(wall_time, datetime.timedelta)
for wall_time in self.data.metadata["wall_time"])
if "cpu_time" in self.data.metadata:
assert self.data.metadata["cpu_time"]
assert all(isinstance(cpu_time, datetime.timedelta)
for cpu_time in self.data.metadata["cpu_time"])
class GenericHFSPTest(GenericSPTest):
# Approximate HF energy of dvb after SCF in STO-3G (from DALTON 2015).
hf_scfenergy = -10334.03948035995
# Approximate energy of the innermost molecular orbital of DVB with
# HF/STO-3G (from Psi4 1.3.1).
hf_moenergy = -300.43401785663235
@skipForParser('FChk', 'Formatted Checkpoint files do not have a section for SCF energy')
def testscfenergy(self):
"""Is the SCF energy within the target?"""
assert abs(self.data.scfenergies[-1]-self.hf_scfenergy) < 6.5e-1
def testfirstmoenergy(self):
"""Is the lowest energy molecular orbital within the target?"""
assert abs(self.data.moenergies[0][0]-self.hf_moenergy) < 1.6e-1
class ADFSPTest(GenericSPTest):
"""Customized restricted single point unittest"""
# ADF only prints up to 0.1mD per atom, so the precision here is worse than 0.1mD.
mass_precision = 0.3
foverlap00 = 1.00003
foverlap11 = 1.02672
foverlap22 = 1.03585
num_scf_criteria = 2
b3lyp_energy = -140
b3lyp_moenergy = -269.6079423873336
def testfoverlaps(self):
"""Are the dims and values of the fragment orbital overlap matrix correct?"""
assert self.data.fooverlaps.shape == (self.data.nbasis, self.data.nbasis)
# The matrix is symmetric.
row = self.data.fooverlaps[0,:]
col = self.data.fooverlaps[:,0]
assert sum(col - row) == 0.0
# Although the diagonal elements are close to zero, the SFOs
# are generally not normalized, so test for a few specific values.
assert abs(self.data.fooverlaps[0, 0]-self.foverlap00) < 0.0001
assert abs(self.data.fooverlaps[1, 1]-self.foverlap11) < 0.0001
assert abs(self.data.fooverlaps[2, 2]-self.foverlap22) < 0.0001
class GaussianSPTest(GenericSPTest):
"""Customized restricted single point unittest"""
num_scf_criteria = 3
class JaguarSPTest(GenericSPTest):
"""Customized restricted single point KS unittest"""
num_scf_criteria = 2
class JaguarHFSPTest(JaguarSPTest, GenericHFSPTest):
"""Customized restricted single point KS unittest"""
class Jaguar7SPTest(JaguarSPTest):
"""Customized restricted single point unittest"""
# Jaguar prints only 10 virtual MOs by default. Can we re-run with full output?
def testlengthmoenergies(self):
"""Is the number of evalues equal to the number of occ. MOs + 10?"""
assert len(self.data.moenergies[0]) == self.data.homos[0]+11
class MolcasSPTest(GenericSPTest):
"""Customized restricted single point unittest"""
num_scf_criteria = 4
class MolproSPTest(GenericHFSPTest):
"""Customized restricted single point HF unittest"""
num_scf_criteria = 2
class NWChemKSSPTest(GenericSPTest):
"""Customized restricted single point unittest"""
num_scf_criteria = 3
class PsiSPTest(GenericSPTest):
"""Customized restricted single point KS unittest"""
num_scf_criteria = 2
class PsiHFSPTest(PsiSPTest, GenericHFSPTest):
"""Customized restricted single point HF unittest"""
class OrcaSPTest(GenericSPTest):
"""Customized restricted single point unittest"""
# Orca has different weights for the masses
molecularmass = 130190
b3lyp_moenergy_delta = 1.2e-1
num_scf_criteria = 3
class TurbomoleSPTest(GenericSPTest):
"""Customized restricted single point KS unittest"""
num_scf_criteria = 2
def testmetadata_basis_set(self):
"""Does metadata have expected keys and values?"""
# One of our test cases used sto-3g hondo
valid_basis = self.data.metadata["basis_set"].lower() in ("sto-3g", "sto-3g hondo")
assert valid_basis
class TurbomoleHFSPTest(TurbomoleSPTest, GenericHFSPTest):
"""Customized restricted single point HF unittest"""
class GenericDispersionTest(unittest.TestCase):
"""Generic single-geometry dispersion correction unittest"""
dispersionenergy = -0.4005496
def testdispersionenergies(self):
"""Is the dispersion energy parsed correctly?"""
assert len(self.data.dispersionenergies) == 1
assert abs(self.data.dispersionenergies[0]-self.dispersionenergy) < 2.0e-7
class FireflyDispersionTest(GenericDispersionTest):
"""Customized single-geometry dispersion correction unittest"""
dispersionenergy = -0.4299821
class SolventMetadataTest(unittest.TestCase):
"""Check we can parse implicit solvent data."""
model = ""
# Toluene
static_dielectric_constant = 2.3741
def test_solvent_model(self) -> None:
"""Check solvent model was parsed correctly"""
assert self.data.metadata['solvent_model'] == self.model
@skipForLogfile("basicQChem6.0/water_hf_solvent_smd_iefpcm.out", "the internally-used dielectric constant isn't printed, only solvent name")
def test_solvent_dielectric(self) -> None:
"""Check solvent dielectric was parsed correctly"""
assert abs(self.data.metadata['solvent_params']['epsilon'] - self.static_dielectric_constant) < 1.0e-4
class QChemSolventMetadataTest(SolventMetadataTest):
static_dielectric_constant = 2.370
class IEFPCMMetadataTest(SolventMetadataTest):
"""Check we can parse implicit solvent data."""
model = "IEFPCM"
class SCIPCMMetadataTest(SolventMetadataTest):
"""Check we can parse implicit solvent data."""
model = "SCIPCM"
class IPCMMetadataTest(SolventMetadataTest):
"""Check we can parse implicit solvent data."""
model = "IPCM"
static_dielectric_constant = 78.3
class COSMOMetadataTest(SolventMetadataTest):
"""Check we can parse implicit solvent data."""
model = "COSMO"
class CPCMMetadataTest(SolventMetadataTest):
"""Check we can parse implicit solvent data."""
model = "CPCM"
class CPCMCOSMOMetadataTest(SolventMetadataTest):
"""Check we can parse implicit solvent data."""
model = "CPCM-COSMO"
class SMDIEFPCMMetadataTest(SolventMetadataTest):
"""Check we can parse implicit solvent data."""
model = "SMD-IEFPCM"
class SMDCPCMMetadataTest(SolventMetadataTest):
"""Check we can parse implicit solvent data."""
model = "SMD-CPCM"
class QChemSMDIEFPCMMetadataTest(QChemSolventMetadataTest, SMDIEFPCMMetadataTest):
"""Check we can parse implicit solvent data."""
class QChemSMDCPCMMetadataTest(QChemSolventMetadataTest, SMDCPCMMetadataTest):
"""Check we can parse implicit solvent data."""
if __name__ == "__main__":
import sys
sys.path.insert(1, os.path.join(__filedir__, ".."))
from test_data import DataSuite
suite = DataSuite(['SP'])
suite.testall()
|