votca 2024.2-dev
Loading...
Searching...
No Matches
regular_grid.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2020 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// Local VOTCA includes
25
26namespace votca {
27namespace xtp {
28
29void Regular_Grid::GridSetup(const Eigen::Array3d& stepsizes,
30 const Eigen::Array3d& padding,
31 const QMMolecule& atoms, const AOBasis& basis) {
32
33 std::pair<Eigen::Vector3d, Eigen::Vector3d> extension =
34 atoms.CalcSpatialMinMax();
35 Eigen::Array3d min = extension.first.array();
36 Eigen::Array3d max = extension.second.array();
37 Eigen::Array3d doublesteps = (max - min + 2 * padding) / stepsizes + 1.0;
38 Eigen::Array<votca::Index, 3, 1> steps = (doublesteps.ceil()).cast<Index>();
39
40 // needed to symmetrize grid around molecule
41 Eigen::Array3d padding_sym =
42 (doublesteps - steps.cast<double>()) * stepsizes * 0.5 + padding;
43 GridSetup(steps, padding_sym, atoms, basis);
44}
45
46void Regular_Grid::GridSetup(const Eigen::Array<Index, 3, 1>& steps,
47 const Eigen::Array3d& padding,
48 const QMMolecule& atoms, const AOBasis& basis) {
49
50 std::pair<Eigen::Vector3d, Eigen::Vector3d> extension =
51 atoms.CalcSpatialMinMax();
52 Eigen::Array3d min = extension.first.array();
53 Eigen::Array3d max = extension.second.array();
54 startingpoint_ = min - padding;
55 stepsizes_ = (max - min + 2 * padding) / (steps - 1).cast<double>();
56 steps_ = steps;
57 const Index gridboxsize = 500;
58
59 GridBox gridbox;
60 for (Index i = 0; i < steps.x(); i++) {
61 double x = startingpoint_.x() + double(i) * stepsizes_.x();
62 for (Index j = 0; j < steps.y(); j++) {
63 double y = startingpoint_.y() + double(j) * stepsizes_.y();
64 for (Index k = 0; k < steps.z(); k++) {
65 double z = startingpoint_.z() + double(k) * stepsizes_.z();
67 point.grid_weight = 1.0;
68 point.grid_pos = Eigen::Vector3d(x, y, z);
69 gridbox.addGridPoint(point);
70 if (gridbox.size() == gridboxsize) {
71 grid_boxes_.push_back(gridbox);
72 gridbox = GridBox();
73 }
74 }
75 }
76 }
77 if (gridbox.size() > 0) {
78 grid_boxes_.push_back(gridbox);
79 }
80#pragma omp parallel for
81 for (Index i = 0; i < getBoxesSize(); i++) {
82 grid_boxes_[i].FindSignificantShells(basis);
83 grid_boxes_[i].PrepareForIntegration();
84 }
85
87 for (auto& box : grid_boxes_) {
88 totalgridsize_ += box.size();
89 }
90}
91
92} // namespace xtp
93} // namespace votca
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
std::pair< Eigen::Vector3d, Eigen::Vector3d > CalcSpatialMinMax() const
Index size() const
Definition gridbox.h:51
void addGridPoint(const GridContainers::Cartesian_gridpoint &point)
Definition gridbox.h:63
void GridSetup(const Eigen::Array< Index, 3, 1 > &steps, const Eigen::Array3d &padding, const QMMolecule &atoms, const AOBasis &basis)
Eigen::Array3d stepsizes_
std::vector< GridBox > grid_boxes_
Eigen::Vector3d startingpoint_
Eigen::Array< Index, 3, 1 > steps_
Index getBoxesSize() const
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26