source: src/Actions/FillAction/SuspendInMoleculeAction.cpp@ 3c9ac3

Action_Thermostats Adding_MD_integration_tests Adding_StructOpt_integration_tests AutomationFragmentation_failures Candidate_v1.6.1 ChemicalSpaceEvaluator Enhanced_StructuralOptimization Enhanced_StructuralOptimization_continued Exclude_Hydrogens_annealWithBondGraph Fix_Verbose_Codepatterns ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_continued ForceAnnealing_with_BondGraph_continued_betteresults ForceAnnealing_with_BondGraph_contraction-expansion Gui_displays_atomic_force_velocity JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool PythonUI_with_named_parameters Recreated_GuiChecks StoppableMakroAction TremoloParser_IncreasedPrecision
Last change on this file since 3c9ac3 was 9eb71b3, checked in by Frederik Heber <frederik.heber@…>, 8 years ago

Commented out MemDebug include and Memory::ignore.

  • MemDebug clashes with various allocation operators that use a specific placement in memory. It is so far not possible to wrap new/delete fully. Hence, we stop this effort which so far has forced us to put ever more includes (with clashes) into MemDebug and thereby bloat compilation time.
  • MemDebug does not add that much usefulness which is not also provided by valgrind.
  • Property mode set to 100644
File size: 12.6 KB
Line 
1/*
2 * Project: MoleCuilder
3 * Description: creates and alters molecular systems
4 * Copyright (C) 2014 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 * SuspendInMoleculeAction.cpp
25 *
26 * Created on: Sep 05, 2014
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 "Actions/UndoRedoHelpers.hpp"
38#include "Atom/atom.hpp"
39#include "Atom/AtomicInfo.hpp"
40#include "Atom/CopyAtoms/CopyAtoms_withBonds.hpp"
41#include "Bond/BondInfo.hpp"
42#include "CodePatterns/Log.hpp"
43#include "Descriptors/MoleculeOrderDescriptor.hpp"
44#include "Element/element.hpp"
45#include "Filling/Cluster.hpp"
46#include "Filling/Filler.hpp"
47#include "Filling/Preparators/BoxFillerPreparator.hpp"
48#include "LinkedCell/linkedcell.hpp"
49#include "LinkedCell/PointCloudAdaptor.hpp"
50#include "molecule.hpp"
51#include "Parser/FormatParserInterface.hpp"
52#include "Parser/FormatParserStorage.hpp"
53#include "Tesselation/boundary.hpp"
54#include "Tesselation/tesselation.hpp"
55#include "World.hpp"
56
57#include <algorithm>
58#include<gsl/gsl_poly.h>
59#include <iostream>
60#include <string>
61#include <vector>
62
63#include "Actions/FillAction/SuspendInMoleculeAction.hpp"
64
65using namespace MoleCuilder;
66
67static double calculateMass(const molecule &_mol)
68{
69 // sum up the atomic masses
70 const double mass = _mol.getAtomSet().totalMass();
71 LOG(2, "DEBUG: Molecule "+_mol.getName()+"'s summed mass is "
72 << setprecision(10) << mass << " atomicmassunit.");
73 return mass;
74}
75
76static double calculateEnvelopeVolume(
77 molecule &_mol,
78 std::vector<double> &_diameters)
79{
80 const bool IsAngstroem = true;
81 class Tesselation *TesselStruct = NULL;
82
83 Boundaries *BoundaryPoints = GetBoundaryPoints(&_mol, TesselStruct);
84 const double * diameters =
85 GetDiametersOfCluster(BoundaryPoints, &_mol, TesselStruct, IsAngstroem);
86 std::copy(&diameters[0], &diameters[3], _diameters.begin());
87 delete[] diameters;
88 PointCloudAdaptor< molecule > cloud(&_mol, _mol.getName());
89 LinkedCell_deprecated *LCList = new LinkedCell_deprecated(cloud, 10.);
90 FindConvexBorder(&_mol, BoundaryPoints, TesselStruct, (const LinkedCell_deprecated *&)LCList, NULL);
91 delete (LCList);
92 delete[] BoundaryPoints;
93
94 // some preparations beforehand
95 const double volume = TesselStruct->getVolumeOfConvexEnvelope(IsAngstroem);
96
97 delete TesselStruct;
98
99 LOG(2, "DEBUG: Molecule "+_mol.getName()+"'s volume is "
100 << setprecision(10) << volume << " angstrom^3.");
101
102 return volume;
103}
104
105// and construct the stuff
106#include "SuspendInMoleculeAction.def"
107#include "Action_impl_pre.hpp"
108/** =========== define the function ====================== */
109ActionState::ptr FillSuspendInMoleculeAction::performCall() {
110 // get the filler molecule
111 std::vector<AtomicInfo> movedatoms;
112 molecule *filler = NULL;
113 {
114 const std::vector< molecule *> molecules = World::getInstance().getSelectedMolecules();
115 if (molecules.size() != 1) {
116 STATUS("No exactly one molecule selected, aborting,");
117 return Action::failure;
118 }
119 filler = *(molecules.begin());
120 }
121 for(molecule::const_iterator iter = const_cast<const molecule *>(filler)->begin();
122 iter != const_cast<const molecule *>(filler)->end(); ++iter)
123 movedatoms.push_back( AtomicInfo(*(*iter)) );
124 LOG(1, "INFO: Chosen molecule has " << filler->size() << " atoms.");
125
126 // center filler's tip at origin
127 filler->CenterEdge();
128
129 std::vector<molecule *> molecules = World::getInstance().getAllMolecules();
130 if (molecules.size() < 2) {
131 STATUS("There must be at least two molecules: filler and to be suspended.");
132 return Action::failure;
133 }
134
135 /// first we need to calculate some volumes and masses
136 double totalmass = 0.;
137 const bool IsAngstroem = true;
138 Vector BoxLengths;
139 double clustervolume = 0.;
140 std::vector<double> GreatestDiameter(NDIM, 0.);
141 for (std::vector<molecule *>::const_iterator iter = molecules.begin();
142 iter != molecules.end(); ++iter)
143 {
144 // skip the filler
145 if (*iter == filler)
146 continue;
147 molecule & mol = **iter;
148 const double mass = calculateMass(mol);
149 totalmass += mass;
150 std::vector<double> diameters(NDIM, 0.);
151 const double volume = calculateEnvelopeVolume(mol, diameters);
152 clustervolume += volume;
153 for (size_t i=0;i<NDIM;++i)
154 GreatestDiameter[i] = std::max(GreatestDiameter[i], diameters[i]);
155 }
156 LOG(1, "INFO: The summed mass is " << setprecision(10)
157 << totalmass << " atomicmassunit.");
158 LOG(1, "INFO: The average density is " << setprecision(10)
159 << totalmass / clustervolume << " atomicmassunit/"
160 << (IsAngstroem ? " angstrom" : " atomiclength") << "^3.");
161 if ( ((totalmass / clustervolume < 1.) && (params.density.get() > 1.))
162 || ((totalmass / clustervolume > 1.) && (params.density.get() < 1.))) {
163 STATUS("Desired and present molecular densities must both be either in [0,1) or in (1, inf).");
164 return Action::failure;
165 }
166
167 // calculate maximum solvent density
168 std::vector<double> fillerdiameters(NDIM, 0.);
169 const double fillervolume = calculateEnvelopeVolume(*filler, fillerdiameters);
170 const double fillermass = calculateMass(*filler);
171 LOG(1, "INFO: The filler's mass is " << setprecision(10)
172 << fillermass << " atomicmassunit, and it's volume is "
173 << fillervolume << (IsAngstroem ? " angstrom" : " atomiclength") << "^3.");
174// const double solventdensity = fillermass / fillervolume;
175
176 /// solve cubic polynomial
177 double cellvolume = 0.;
178 LOG(1, "Solving equidistant suspension in water problem ...");
179 // s = solvent, f = filler, 0 = initial molecules/cluster
180 // v_s = v_0 + v_f, m_s = m_0 + rho_f * v_f --> rho_s = m_s/v_s ==> v_f = (m_0 - rho_s * v_o) / (rho_s - rho_f)
181 cellvolume = (totalmass - params.density.get() * clustervolume) / (params.density.get() - 1.) + clustervolume;
182 LOG(1, "Cellvolume needed for a density of " << params.density.get()
183 << " g/cm^3 is " << cellvolume << " angstroem^3.");
184
185 const double minimumvolume =
186 (GreatestDiameter[0] * GreatestDiameter[1] * GreatestDiameter[2]);
187 LOG(1, "Minimum volume of the convex envelope contained in a rectangular box is "
188 << minimumvolume << " angstrom^3.");
189
190 if (minimumvolume > cellvolume) {
191 ELOG(1, "The containing box already has a greater volume than the envisaged cell volume!");
192 LOG(0, "Setting Box dimensions to minimum possible, the greatest diameters.");
193 for (int i = 0; i < NDIM; i++)
194 BoxLengths[i] = GreatestDiameter[i];
195// mol->CenterEdge();
196 } else {
197 BoxLengths[0] = GreatestDiameter[0] + GreatestDiameter[1] + GreatestDiameter[2];
198 BoxLengths[1] = GreatestDiameter[0] * GreatestDiameter[1]
199 + GreatestDiameter[0] * GreatestDiameter[2]
200 + GreatestDiameter[1] * GreatestDiameter[2];
201 BoxLengths[2] = minimumvolume - cellvolume;
202 std::vector<double> x(3, 0.);
203 // for cubic polynomial there are either 1 or 3 unique solutions
204 if (gsl_poly_solve_cubic(BoxLengths[0], BoxLengths[1], BoxLengths[2], &x[0], &x[1], &x[2]) == 1) {
205 x[1] = x[0];
206 x[2] = x[0];
207 } else {
208 std::swap(x[0], x[2]); // sorted in ascending order
209 }
210 LOG(0, "RESULT: The resulting spacing is: " << x << " .");
211
212 cellvolume = 1.;
213 for (size_t i = 0; i < NDIM; ++i) {
214 BoxLengths[i] = x[i] + GreatestDiameter[i];
215 cellvolume *= BoxLengths[i];
216 }
217 }
218
219 // TODO: Determine counts from resulting mass correctly (hard problem due to integers)
220 std::vector<unsigned int> counts(3, 0);
221 const unsigned int totalcounts = round(params.density.get() * cellvolume - totalmass) / fillermass;
222 if (totalcounts > 0) {
223 counts[0] = ceil(BoxLengths[0]/3.1);
224 counts[1] = ceil(BoxLengths[1]/3.1);
225 counts[2] = ceil(BoxLengths[2]/3.1);
226 }
227
228 // update Box of atoms by boundary
229 {
230 RealSpaceMatrix domain;
231 for(size_t i =0; i<NDIM;++i)
232 domain.at(i,i) = BoxLengths[i];
233 World::getInstance().setDomain(domain);
234 }
235 LOG(0, "RESULT: The resulting cell dimensions are: " << BoxLengths[0] << " and " << BoxLengths[1] << " and " << BoxLengths[2] << " with total volume of " << cellvolume << " " << (IsAngstroem ? "angstrom" : "atomiclength") << "^3.");
236
237 // prepare the filler preparator
238 BoxFillerPreparator filler_preparator(filler);
239 filler_preparator.addVoidPredicate(params.mindistance.get());
240 filler_preparator.addRandomInserter(
241 params.RandAtomDisplacement.get(),
242 params.RandMoleculeDisplacement.get(),
243 params.DoRotate.get());
244 Vector offset(.5,.5,.5);
245 filler_preparator.addCubeMesh(
246 counts,
247 offset,
248 World::getInstance().getDomain().getM());
249 if (!filler_preparator()) {
250 STATUS("Filler was not fully constructed.");
251 return Action::failure;
252 }
253
254 // use filler
255 bool successflag = false;
256 FillSuspendInMoleculeState *UndoState = NULL;
257 {
258 // fill
259 Filler *fillerFunction = filler_preparator.obtainFiller();
260 // TODO: When molecule::getBoundingSphere() does not use a sphere anymore,
261 // we need to check whether we rotate the molecule randomly. For this to
262 // work we need a sphere!
263 const Shape s = filler->getBoundingSphere(params.RandAtomDisplacement.get());
264 ClusterInterface::Cluster_impl cluster( new Cluster(filler->getAtomIds(), s) );
265 CopyAtoms_withBonds copyMethod;
266 Filler::ClusterVector_t ClonedClusters;
267 successflag = (*fillerFunction)(copyMethod, cluster, ClonedClusters);
268 delete fillerFunction;
269
270 // append each cluster's atoms to clonedatoms (however not selected ones)
271 std::vector<const atom *> clonedatoms;
272 std::vector<AtomicInfo> clonedatominfos;
273 for (Filler::ClusterVector_t::const_iterator iter = ClonedClusters.begin();
274 iter != ClonedClusters.end(); ++iter) {
275 const AtomIdSet &atoms = (*iter)->getAtomIds();
276 clonedatoms.reserve(clonedatoms.size()+atoms.size());
277 for (AtomIdSet::const_iterator atomiter = atoms.begin(); atomiter != atoms.end(); ++atomiter)
278 if (!filler->containsAtom(*atomiter)) {
279 clonedatoms.push_back( *atomiter );
280 clonedatominfos.push_back( AtomicInfo(*(*atomiter)) );
281 }
282 }
283 std::vector< BondInfo > clonedbonds;
284 StoreBondInformationFromAtoms(clonedatoms, clonedbonds);
285 LOG(2, "DEBUG: There are " << clonedatominfos.size() << " newly created atoms.");
286
287 if (!successflag) {
288 STATUS("Insertion failed, removing inserted clusters, translating original one back");
289 RemoveAtomsFromAtomicInfo(clonedatominfos);
290 clonedatoms.clear();
291 SetAtomsFromAtomicInfo(movedatoms);
292 } else {
293 std::vector<Vector> MovedToVector(filler->size(), zeroVec);
294 std::transform(filler->begin(), filler->end(), MovedToVector.begin(),
295 boost::bind(&AtomInfo::getPosition, _1) );
296 UndoState = new FillSuspendInMoleculeState(clonedatominfos,clonedbonds,movedatoms,MovedToVector,params);
297 }
298 }
299
300 if (successflag)
301 return ActionState::ptr(UndoState);
302 else {
303 return Action::failure;
304 }
305}
306
307ActionState::ptr FillSuspendInMoleculeAction::performUndo(ActionState::ptr _state) {
308 FillSuspendInMoleculeState *state = assert_cast<FillSuspendInMoleculeState*>(_state.get());
309
310 // remove all created atoms
311 RemoveAtomsFromAtomicInfo(state->clonedatoms);
312 // add the original cluster
313 SetAtomsFromAtomicInfo(state->movedatoms);
314
315 return ActionState::ptr(_state);
316}
317
318ActionState::ptr FillSuspendInMoleculeAction::performRedo(ActionState::ptr _state){
319 FillSuspendInMoleculeState *state = assert_cast<FillSuspendInMoleculeState*>(_state.get());
320
321 // place filler cluster again at new spot
322 ResetAtomPosition(state->movedatoms, state->MovedToVector);
323
324 // re-create all clusters
325 bool statusflag = AddAtomsFromAtomicInfo(state->clonedatoms);
326
327 // re-create the bonds
328 if (statusflag)
329 AddBondsFromBondInfo(state->clonedbonds);
330 if (statusflag)
331 return ActionState::ptr(_state);
332 else {
333 STATUS("Failed re-adding filled in atoms.");
334 return Action::failure;
335 }
336}
337
338bool FillSuspendInMoleculeAction::canUndo() {
339 return false;
340}
341
342bool FillSuspendInMoleculeAction::shouldUndo() {
343 return false;
344}
345/** =========== end of function ====================== */
Note: See TracBrowser for help on using the repository browser.