votca 2024.1-dev
Loading...
Searching...
No Matches
ecpaobasis.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2021 The VOTCA Development Team
3 * (http://www.votca.org)
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License")
6 *
7 * You may not use this file except in compliance with the License.
8 * You may obtain a copy of the License at
9 *
10 * http://www.apache.org/licenses/LICENSE-2.0
11 *
12 * Unless required by applicable law or agreed to in writing, software
13 * distributed under the License is distributed on an "AS IS" BASIS,
14 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
15 * See the License for the specific language governing permissions and
16 * limitations under the License.
17 *
18 */
19
20#pragma once
21#ifndef VOTCA_XTP_ECPAOBASIS_H
22#define VOTCA_XTP_ECPAOBASIS_H
23// Local VOTCA includes
24#include "checkpoint.h"
25#include "ecpbasisset.h"
26#include "qmatom.h"
27#include <libecpint/ecp.hpp>
28
29namespace votca {
30namespace xtp {
31class QMMolecule;
32class ECPBasisSet;
33class CheckpointWriter;
34class CheckpointReader;
35
36std::ostream& operator<<(std::ostream& out, const libecpint::ECP& ecp);
37
44 public:
45 // returns element names for which no ecp was found
46 std::vector<std::string> Fill(const ECPBasisSet& bs, QMMolecule& atoms);
47
48 using constECPAOShellIterator = std::vector<libecpint::ECP>::const_iterator;
49 constECPAOShellIterator begin() const { return aopotentials_.begin(); }
50 constECPAOShellIterator end() const { return aopotentials_.end(); }
51
52 using ECPAOShellIterator = std::vector<libecpint::ECP>::iterator;
55
56 Index getMaxL() const;
57 void AddECPChargeToMolecule(QMMolecule& mol) const;
58
59 const std::string& Name() const { return name_; }
60
61 void UpdatePotentialPositions(const QMMolecule& mol);
62
63 void WriteToCpt(CheckpointWriter& w) const;
64
66
67 void add(const ECPAOBasis& other);
68
69 friend std::ostream& operator<<(std::ostream& out, const ECPAOBasis& ecp);
70
71 private:
72 void clear();
73
74 std::vector<Index> ncore_perAtom_;
75
76 std::vector<libecpint::ECP> aopotentials_;
77
78 std::string name_ = "";
79};
80
81} // namespace xtp
82} // namespace votca
83
84#endif // VOTCA_XTP_ECPAOBASIS_H
Container to hold ECPs for all atoms.
Definition ecpaobasis.h:43
void ReadFromCpt(CheckpointReader &r)
constECPAOShellIterator begin() const
Definition ecpaobasis.h:49
Index getMaxL() const
Definition ecpaobasis.cc:31
std::vector< libecpint::ECP > aopotentials_
Definition ecpaobasis.h:76
std::vector< libecpint::ECP >::iterator ECPAOShellIterator
Definition ecpaobasis.h:52
constECPAOShellIterator end() const
Definition ecpaobasis.h:50
friend std::ostream & operator<<(std::ostream &out, const ECPAOBasis &ecp)
const std::string & Name() const
Definition ecpaobasis.h:59
void UpdatePotentialPositions(const QMMolecule &mol)
Definition ecpaobasis.cc:51
std::vector< std::string > Fill(const ECPBasisSet &bs, QMMolecule &atoms)
Definition ecpaobasis.cc:69
void WriteToCpt(CheckpointWriter &w) const
ECPAOShellIterator begin()
Definition ecpaobasis.h:53
std::vector< Index > ncore_perAtom_
Definition ecpaobasis.h:74
void AddECPChargeToMolecule(QMMolecule &mol) const
Definition ecpaobasis.cc:39
void add(const ECPAOBasis &other)
Definition ecpaobasis.cc:58
ECPAOShellIterator end()
Definition ecpaobasis.h:54
std::vector< libecpint::ECP >::const_iterator constECPAOShellIterator
Definition ecpaobasis.h:48
std::ostream & operator<<(std::ostream &out, const Correlate &c)
Definition correlate.h:53
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26