source: src/Fragmentation/Interfragmenter.cpp@ d45ed9

Action_Thermostats Add_AtomRandomPerturbation Add_FitFragmentPartialChargesAction Add_RotateAroundBondAction Add_SelectAtomByNameAction Added_ParseSaveFragmentResults Adding_Graph_to_ChangeBondActions Adding_MD_integration_tests Adding_ParticleName_to_Atom Adding_StructOpt_integration_tests AtomFragments Automaking_mpqc_open AutomationFragmentation_failures Candidate_v1.5.4 Candidate_v1.6.0 Candidate_v1.6.1 ChangeBugEmailaddress ChangingTestPorts ChemicalSpaceEvaluator CombiningParticlePotentialParsing Combining_Subpackages Debian_Package_split Debian_package_split_molecuildergui_only Disabling_MemDebug Docu_Python_wait EmpiricalPotential_contain_HomologyGraph EmpiricalPotential_contain_HomologyGraph_documentation Enable_parallel_make_install Enhance_userguide Enhanced_StructuralOptimization Enhanced_StructuralOptimization_continued Example_ManyWaysToTranslateAtom Exclude_Hydrogens_annealWithBondGraph FitPartialCharges_GlobalError Fix_ChargeSampling_PBC Fix_ChronosMutex Fix_FitPartialCharges Fix_FitPotential_needs_atomicnumbers Fix_ForceAnnealing Fix_IndependentFragmentGrids Fix_ParseParticles Fix_ParseParticles_split_forward_backward_Actions Fix_PopActions Fix_QtFragmentList_sorted_selection Fix_Restrictedkeyset_FragmentMolecule Fix_StatusMsg Fix_StepWorldTime_single_argument Fix_Verbose_Codepatterns Fix_fitting_potentials Fixes ForceAnnealing_goodresults ForceAnnealing_oldresults ForceAnnealing_tocheck ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_continued ForceAnnealing_with_BondGraph_continued_betteresults ForceAnnealing_with_BondGraph_contraction-expansion FragmentAction_writes_AtomFragments FragmentMolecule_checks_bonddegrees GeometryObjects Gui_Fixes Gui_displays_atomic_force_velocity IndependentFragmentGrids IndependentFragmentGrids_IndividualZeroInstances IndependentFragmentGrids_IntegrationTest IndependentFragmentGrids_Sole_NN_Calculation JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool JobMarket_unresolvable_hostname_fix ODR_violation_mpqc_open PartialCharges_OrthogonalSummation PdbParser_setsAtomName PythonUI_with_named_parameters QtGui_reactivate_TimeChanged_changes Recreated_GuiChecks Rewrite_FitPartialCharges RotateToPrincipalAxisSystem_UndoRedo SaturateAtoms_findBestMatching SaturateAtoms_singleDegree StoppableMakroAction Subpackage_CodePatterns Subpackage_JobMarket Subpackage_LinearAlgebra Subpackage_levmar Subpackage_mpqc_open Subpackage_vmg Switchable_LogView ThirdParty_MPQC_rebuilt_buildsystem TrajectoryDependenant_MaxOrder TremoloParser_IncreasedPrecision TremoloParser_MultipleTimesteps TremoloParser_setsAtomName Ubuntu_1604_changes stable
Last change on this file since d45ed9 was d45ed9, checked in by Frederik Heber <heber@…>, 9 years ago

Abstracted atom -> fragment association out of Interfragmenter.

  • Property mode set to 100644
File size: 9.2 KB
Line 
1/*
2 * Project: MoleCuilder
3 * Description: creates and alters molecular systems
4 * Copyright (C) 2013 Frederik Heber. All rights reserved.
5 *
6 *
7 * This file is part of MoleCuilder.
8 *
9 * MoleCuilder is free software: you can redistribute it and/or modify
10 * it under the terms of the GNU General Public License as published by
11 * the Free Software Foundation, either version 2 of the License, or
12 * (at your option) any later version.
13 *
14 * MoleCuilder is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 * GNU General Public License for more details.
18 *
19 * You should have received a copy of the GNU General Public License
20 * along with MoleCuilder. If not, see <http://www.gnu.org/licenses/>.
21 */
22
23/*
24 * Interfragmenter.cpp
25 *
26 * Created on: Jul 5, 2013
27 * Author: heber
28 */
29
30// include config.h
31#ifdef HAVE_CONFIG_H
32#include <config.h>
33#endif
34
35#include "CodePatterns/MemDebug.hpp"
36
37#include "Interfragmenter.hpp"
38
39#include "CodePatterns/Assert.hpp"
40#include "CodePatterns/Log.hpp"
41
42#include "LinearAlgebra/Vector.hpp"
43
44#include "Descriptors/AtomDescriptor.hpp"
45#include "Element/element.hpp"
46#include "Fragmentation/Graph.hpp"
47#include "Fragmentation/KeySet.hpp"
48#include "LinkedCell/LinkedCell_View.hpp"
49#include "LinkedCell/types.hpp"
50#include "World.hpp"
51
52static Vector getAtomIdSetCenter(
53 const AtomIdSet &_atoms)
54{
55 const molecule * const _mol = (*_atoms.begin())->getMolecule();
56 const size_t atoms_size = _atoms.getAtomIds().size();
57 Vector center;
58 for (AtomIdSet::const_iterator iter = _atoms.begin();
59 iter != _atoms.end(); ++iter) {
60 center += (*iter)->getPosition();
61 ASSERT ( _mol == (*iter)->getMolecule(),
62 "getAtomIdSetCenter() - ids in same keyset belong to different molecule.");
63 }
64 center *= 1./(double)atoms_size;
65
66 return center;
67}
68
69Interfragmenter::candidates_t Interfragmenter::getNeighborsOutsideMolecule(
70 const AtomIdSet &_atoms,
71 double _Rcut,
72 const enum HydrogenTreatment _treatment) const
73{
74 /// go through linked cell and get all neighboring atoms up to Rcut
75 const LinkedCell::LinkedCell_View view = World::getInstance().getLinkedCell(_Rcut);
76 const Vector center = getAtomIdSetCenter(_atoms);
77 const LinkedCell::LinkedList neighbors = view.getAllNeighbors(_Rcut, center);
78 LOG(4, "DEBUG: Obtained " << neighbors.size() << " neighbors in distance of "
79 << _Rcut << " around " << center << ".");
80
81 /// remove all atoms that belong to same molecule as does one of the
82 /// fragment's atoms
83 const molecule * const _mol = (*_atoms.begin())->getMolecule();
84 candidates_t candidates;
85 candidates.reserve(neighbors.size());
86 for (LinkedCell::LinkedList::const_iterator iter = neighbors.begin();
87 iter != neighbors.end(); ++iter) {
88 const atom * const _atom = static_cast<const atom * const >(*iter);
89 ASSERT( _atom != NULL,
90 "Interfragmenter::getNeighborsOutsideMolecule() - a neighbor is not actually an atom?");
91 if ((_atom->getMolecule() != _mol)
92 && (_atom->getPosition().DistanceSquared(center) < _Rcut*_Rcut)
93 && ((_treatment == IncludeHydrogen) || (_atom->getType()->getAtomicNumber() != 1))) {
94 candidates.push_back(_atom);
95 }
96 }
97 LOG(3, "DEBUG: There remain " << candidates.size() << " candidates.");
98
99 return candidates;
100}
101
102Interfragmenter::atomkeyset_t Interfragmenter::getCandidatesSpecificKeySetMap(
103 const candidates_t &_candidates,
104 const atomkeyset_t &_atomkeyset) const
105{
106 atomkeyset_t fragmentmap;
107 for (candidates_t::const_iterator candidateiter = _candidates.begin();
108 candidateiter != _candidates.end(); ++candidateiter) {
109 const atom * _atom = *candidateiter;
110 atomkeyset_t::const_iterator iter = _atomkeyset.find(_atom);
111 ASSERT( iter != _atomkeyset.end(),
112 "Interfragmenter::getAtomSpecificKeySetMap() - could not find atom "
113 +toString(_atom->getId())+" in lookup.");
114 fragmentmap.insert( std::make_pair( _atom, iter->second) );
115 }
116 LOG(4, "DEBUG: Copied part of lookup map contains " << fragmentmap.size() << " keys.");
117
118 return fragmentmap;
119}
120
121void Interfragmenter::combineFragments(
122 const size_t _MaxOrder,
123 const candidates_t &_candidates,
124 atomkeyset_t &_fragmentmap,
125 const KeySet &_keyset,
126 Graph &_InterFragments,
127 int &_counter)
128{
129 for (candidates_t::const_iterator candidateiter = _candidates.begin();
130 candidateiter != _candidates.end(); ++candidateiter) {
131 const atom *_atom = *candidateiter;
132 LOG(3, "DEBUG: Current candidate is " << *_atom << ".");
133 atomkeyset_t::iterator finditer = _fragmentmap.find(_atom);
134 ASSERT( finditer != _fragmentmap.end(),
135 "Interfragmenter::combineFragments() - could not find atom "
136 +toString(_atom->getId())+" in fragment specific lookup.");
137 keysets_t &othersets = finditer->second;
138 ASSERT( !othersets.empty(),
139 "Interfragmenter::combineFragments() - keysets to "+toString(_atom->getId())+
140 "is empty.");
141 keysets_t::iterator otheriter = othersets.begin();
142 while (otheriter != othersets.end()) {
143 const KeySet &otherset = **otheriter;
144 LOG(3, "DEBUG: Current keyset is " << otherset << ".");
145 // only add them one way round and not the other
146 if (otherset < _keyset) {
147 ++otheriter;
148 continue;
149 }
150 // only add if combined they don't exceed the desired maxorder
151 if (otherset.size() + _keyset.size() > _MaxOrder) {
152 LOG(3, "INFO: Rejecting " << otherset << " as in sum their orders exceed " << _MaxOrder);
153 ++otheriter;
154 continue;
155 }
156 KeySet newset(otherset);
157 newset.insert(_keyset.begin(), _keyset.end());
158 LOG(3, "DEBUG: Inserting new combined set " << newset << ".");
159 _InterFragments.insert( std::make_pair(newset, std::make_pair(++_counter, 1.)));
160 // finally, remove the set such that no other combination exists
161 otheriter = othersets.erase(otheriter);
162 }
163 }
164}
165
166void Interfragmenter::operator()(
167 const size_t MaxOrder,
168 const double Rcut,
169 const enum HydrogenTreatment treatment)
170{
171 AtomFragmentsMap AtomFragments(TotalGraph, MaxOrder);
172 const atomkeyset_t &atomkeyset = AtomFragments.getMap();
173
174 Graph InterFragments;
175 int counter = TotalGraph.size();
176
177 /// go through all fragments up to MaxOrder
178 LOG(1, "INFO: Creating inter-fragments.");
179 for (Graph::const_iterator keysetiter = TotalGraph.begin();
180 keysetiter != TotalGraph.end(); ++keysetiter) {
181 const KeySet &keyset = keysetiter->first;
182 LOG(2, "DEBUG: Current keyset is " << keyset);
183 const AtomIdSet atoms(keyset);
184 const size_t atoms_size = atoms.getAtomIds().size();
185 if ((atoms_size > MaxOrder) || (atoms_size == 0))
186 continue;
187
188 // get neighboring atoms outside the current molecule
189 candidates_t candidates = getNeighborsOutsideMolecule(atoms, Rcut, treatment);
190
191 // create a lookup specific to this fragment
192 atomkeyset_t fragmentmap = getCandidatesSpecificKeySetMap(candidates, atomkeyset);
193
194 /// combine each remaining fragment with current fragment to a new fragment
195 /// if keyset is less (to prevent addding same inter-fragment twice)
196 combineFragments(MaxOrder, candidates, fragmentmap, keyset, InterFragments, counter);
197 }
198
199 /// eventually, add all new fragments to the Graph
200 counter = TotalGraph.size();
201 TotalGraph.InsertGraph(InterFragments, counter);
202}
203
204double Interfragmenter::findLargestCutoff(
205 const size_t _MaxOrder,
206 const double _upperbound,
207 const enum HydrogenTreatment _treatment) const
208{
209 double Rcut = _upperbound*_upperbound;
210 std::pair<atomId_t, atomId_t> ClosestPair;
211 // place all atoms into LC grid with some upper bound
212 const LinkedCell::LinkedCell_View view = World::getInstance().getLinkedCell(_upperbound);
213
214 // go through each atom and find closest atom not in the same keyset
215 for (Graph::const_iterator keysetiter = TotalGraph.begin();
216 keysetiter != TotalGraph.end(); ++keysetiter) {
217 const KeySet &keyset = keysetiter->first;
218 LOG(2, "DEBUG: Current keyset is " << keyset);
219 const AtomIdSet atoms(keyset);
220
221 // get neighboring atoms outside the current molecule
222 const candidates_t candidates = getNeighborsOutsideMolecule(atoms, _upperbound, _treatment);
223 const Vector center = getAtomIdSetCenter(atoms);
224
225 for (candidates_t::const_iterator candidateiter = candidates.begin();
226 candidateiter != candidates.end(); ++candidateiter) {
227 atom const * const Walker = *candidateiter;
228 // go through each atom in set and pick minimal distance
229 for (AtomIdSet::const_iterator setiter = atoms.begin(); setiter != atoms.end(); ++setiter) {
230 const double distanceSquared = Walker->getPosition().DistanceSquared(center);
231 // pick the smallest compared to current Rcut if smaller
232 if (distanceSquared < Rcut) {
233 Rcut = distanceSquared;
234 ClosestPair.first = (*setiter)->getId();
235 ClosestPair.second = Walker->getId();
236 LOG(2, "DEBUG: Found new pair " << ClosestPair << " with distance " << sqrt(Rcut));
237 }
238 }
239 }
240 }
241 const double largest_distance = sqrt(Rcut);
242 LOG(1, "INFO: Largest inter-fragment distance to cause no additional fragments: "
243 << largest_distance);
244
245 return largest_distance;
246}
Note: See TracBrowser for help on using the repository browser.