/*
* Project: MoleCuilder
* Description: creates and alters molecular systems
* Copyright (C) 2016 Frederik Heber. All rights reserved.
*
*
* This file is part of MoleCuilder.
*
* MoleCuilder is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 2 of the License, or
* (at your option) any later version.
*
* MoleCuilder is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with MoleCuilder. If not, see .
*/
/*
* FitFragmentPartialChargesAction.cpp
*
* Created on: Oct 09, 2016
* Author: heber
*/
// include config.h
#ifdef HAVE_CONFIG_H
#include
#endif
// needs to come before MemDebug due to placement new
#include
#include "CodePatterns/MemDebug.hpp"
#include "Atom/atom.hpp"
#include "CodePatterns/Log.hpp"
#include "Fragmentation/Exporters/ExportGraph_ToFiles.hpp"
#include "Fragmentation/Graph.hpp"
#include "World.hpp"
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include "Actions/PotentialAction/FitFragmentPartialChargesAction.hpp"
#include "Potentials/PartialNucleiChargeFitter.hpp"
#include "AtomIdSet.hpp"
#include "Descriptors/AtomIdDescriptor.hpp"
#include "Element/element.hpp"
#include "Element/periodentafel.hpp"
#include "Fragmentation/Homology/AtomFragmentsMap.hpp"
#include "Fragmentation/Homology/HomologyGraph.hpp"
#include "Fragmentation/Summation/Containers/FragmentationResultContainer.hpp"
#include "Fragmentation/Summation/Containers/PartialChargesFused.hpp"
#include "Fragmentation/Summation/Containers/PartialChargesMap.hpp"
#include "Fragmentation/Summation/SetValues/IndexedPartialCharges.hpp"
#include "Fragmentation/Summation/SetValues/SamplingGrid.hpp"
#include "Fragmentation/Summation/ZeroInstanceInitializer.hpp"
#include "FunctionApproximation/Extractors.hpp"
#include "Potentials/PartialNucleiChargeFitter.hpp"
#include "Potentials/Particles/ParticleFactory.hpp"
#include "Potentials/Particles/ParticleRegistry.hpp"
#include "Potentials/SerializablePotential.hpp"
#include "World.hpp"
using namespace MoleCuilder;
// and construct the stuff
#include "FitFragmentPartialChargesAction.def"
#include "Action_impl_pre.hpp"
/** =========== define the function ====================== */
namespace detail {
typedef std::map KeysetsToGraph_t;
typedef std::map GraphFittedChargeMap_t;
typedef std::map fitted_charges_t;
typedef std::map GraphIndex_t;
typedef std::set AtomsGraphIndices_t;
typedef boost::bimaps::bimap<
boost::bimaps::multiset_of,
atomId_t > GraphIndices_t;
typedef std::map, std::map > AtomParticleNames_t;
typedef std::map, std::string> GraphToNameMap_t;
typedef std::map fittedcharges_per_fragment_t;
typedef std::map indices_to_keysets_t;
typedef std::map< std::string, std::pair > average_charge_map_t;
};
static void enforceZeroTotalCharge(
PartialNucleiChargeFitter::charges_t &_averaged_charges)
{
double charge_sum = 0.;
charge_sum = std::accumulate(_averaged_charges.begin(), _averaged_charges.end(), charge_sum);
if (fabs(charge_sum) > MYEPSILON) {
std::transform(
_averaged_charges.begin(), _averaged_charges.end(),
_averaged_charges.begin(),
boost::bind(std::minus(), _1, charge_sum/_averaged_charges.size()));
}
charge_sum = 0.;
charge_sum = std::accumulate(_averaged_charges.begin(), _averaged_charges.end(), charge_sum);
ASSERT( fabs(charge_sum) < MYEPSILON,
"enforceZeroTotalCharge() - enforcing neutral net charge failed, "
+toString(charge_sum)+" is the remaining net charge.");
LOG(2, "DEBUG: final charges with net zero charge are " << _averaged_charges);
}
static
std::set accumulateKeySetsForAtoms(
const AtomFragmentsMap::AtomFragmentsMap_t &_atommap,
const std::vector &_selected_atoms)
{
std::set fragments;
for (std::vector::const_iterator iter = _selected_atoms.begin();
iter != _selected_atoms.end(); ++iter) {
const atomId_t atomid = (*iter)->getId();
const AtomFragmentsMap::AtomFragmentsMap_t::const_iterator atomiter =
_atommap.find(atomid);
if ((*iter)->getElementNo() != 1) {
if (atomiter == _atommap.end()) {
ELOG(2, "There are no fragments associated to " << atomid << ".");
continue;
}
const AtomFragmentsMap::keysets_t &keysets = atomiter->second;
LOG(2, "DEBUG: atom " << atomid << " has " << keysets.size() << " fragments.");
fragments.insert( keysets.begin(), keysets.end() );
} else {
LOG(3, "DEBUG: Skipping atom " << atomid << " as it's hydrogen.");
}
}
return fragments;
}
static
void getKeySetsToGraphMapping(
detail::KeysetsToGraph_t &_keyset_graphs,
const std::set &_fragments,
const AtomFragmentsMap &_atomfragments)
{
for (std::set::const_iterator fragmentiter = _fragments.begin();
fragmentiter != _fragments.end(); ++fragmentiter) {
const KeySet &keyset = *fragmentiter;
const AtomFragmentsMap::indices_t &forceindices = _atomfragments.getFullKeyset(keyset);
ASSERT( !forceindices.empty(),
"getKeySetsToGraphMapping() - force keyset to "+toString(keyset)+" is empty.");
KeySet forcekeyset;
forcekeyset.insert(forceindices.begin(), forceindices.end());
forcekeyset.erase(-1);
const HomologyGraph graph(forcekeyset);
LOG(2, "DEBUG: Associating keyset " << forcekeyset << " with graph " << graph);
_keyset_graphs.insert( std::make_pair(keyset, graph) );
}
}
static
bool getPartialChargesForAllIndexSets(
detail::fittedcharges_per_fragment_t &_fittedcharges_per_fragment,
const FragmentationResultContainer &_container,
const std::set &_fragments,
const double _mask_radius,
const bool enforceZeroCharge)
{
const std::map &shortrangedata = _container.getShortRangeResults();
const std::map &longrangedata = _container.getLongRangeResults();
const KeySetsContainer &keysets = _container.getKeySets();
std::map::const_iterator shortiter = shortrangedata.begin();
std::map::const_iterator longiter = longrangedata.begin();
KeySetsContainer::ArrayOfIntVectors::const_iterator keysetiter = keysets.KeySets.begin();
for (; keysetiter != keysets.KeySets.end(); ++keysetiter, ++shortiter, ++longiter) {
// check whether its an allowed keyset (atomfragments has keyset, not forcekeyset)
KeySet currentset;
currentset.insert(keysetiter->begin(), keysetiter->end());
if (_fragments.find(currentset) == _fragments.end())
continue;
// obtain positions and sampled potential
const Fragment::positions_t &fragmentpositions = shortiter->second.positions;
const SamplingGrid &potential = longiter->second.both_sampled_potential;
if ((potential.level == 0)
|| ((potential.begin[0] == potential.end[0])
&& (potential.begin[1] == potential.end[1])
&& (potential.begin[2] == potential.end[2]))) {
ELOG(1, "Sampled grid contains grid made of zero points.");
return 0;
}
// convert std::vector to Vector
PartialNucleiChargeFitter::positions_t positions;
positions.reserve(fragmentpositions.size());
BOOST_FOREACH( Fragment::position_t pos, fragmentpositions) {
positions.push_back( Vector(pos[0], pos[1], pos[2]) );
}
// fit charges
PartialNucleiChargeFitter fitter(potential, positions, _mask_radius);
fitter();
PartialNucleiChargeFitter::charges_t return_charges = fitter.getSolutionAsCharges_t();
LOG(2, "DEBUG: fitted charges are " << return_charges);
// make sum of charges zero if desired
if (enforceZeroCharge)
enforceZeroTotalCharge(return_charges);
// output status info fitted charges
LOG(2, "DEBUG: For fragment " << *keysetiter
<< " we have fitted the following charges " << return_charges << ".");
IndexedPartialCharges::values_t values;
for (PartialNucleiChargeFitter::charges_t::const_iterator chargeiter = return_charges.begin();
chargeiter != return_charges.end(); ++chargeiter) {
values.push_back(partial_charge_t(*chargeiter));
}
_fittedcharges_per_fragment.insert( std::make_pair(currentset, values) );
}
return true;
}
static const atom * getNonHydrogenSurrogate(const atom * const _walker)
{
const atom * surrogate = _walker;
if (surrogate->getElementNo() == 1) {
// it's hydrogen, check its bonding and use its bond partner instead to request
// keysets
const BondList &ListOfBonds = surrogate->getListOfBonds();
if ( ListOfBonds.size() != 1) {
ELOG(1, "Solitary hydrogen in atom " << surrogate->getId() << " detected.");
return _walker;
}
surrogate = (*ListOfBonds.begin())->GetOtherAtom(surrogate);
}
return surrogate;
}
static void addToParticleRegistry(
const ParticleFactory &factory,
ParticleRegistry ®istry,
const periodentafel &periode,
const detail::fitted_charges_t &_fitted_charges,
const detail::GraphIndices_t &_GraphIndices,
detail::AtomParticleNames_t &_atom_particlenames)
{
/// We step here through the atoms and check for all atoms of the
/// same element and which belong to the same set of graphs (via the set
/// of their unique enumeration indices) then they get the same particle
detail::average_charge_map_t average_charge_map;
for (detail::fitted_charges_t::const_iterator chargeiter = _fitted_charges.begin();
chargeiter != _fitted_charges.end(); ++chargeiter) {
const atomId_t &atomid = chargeiter->first;
const double &charge = chargeiter->second;
const atom * const walker = World::getInstance().getAtom(AtomById(atomid));
ASSERT( walker != NULL,
"addToParticleRegistry() - atom "+toString(atomid)
+" not present in the World?");
const detail::GraphIndices_t::right_const_iterator graphiter =
_GraphIndices.right.find(atomid);
ASSERT(graphiter != _GraphIndices.right.end(),
"addToParticleRegistry() - atom #"+toString(atomid)
+" not contained in GraphIndices.");
const detail::AtomParticleNames_t::iterator nameiter =
_atom_particlenames.find(graphiter->second);
const atomicNumber_t elementno = walker->getElementNo();
std::string name;
if ((nameiter != _atom_particlenames.end()) && (nameiter->second.count(elementno))) {
name = (nameiter->second)[elementno];
detail::average_charge_map_t::iterator averageiter = average_charge_map.find(name);
ASSERT( averageiter != average_charge_map.end(),
"addToParticleRegistry() - could not find name "+toString(nameiter->second)
+" in average charge map.");
// add charge to particle
averageiter->second.first += charge;
++(averageiter->second.second);
LOG(1, "INFO: Adding additional charge " << charge << " to particle " << name << " for atom "
<< *walker << " with element " << elementno << ".");
} else {
if (nameiter == _atom_particlenames.end())
_atom_particlenames.insert(
std::make_pair(graphiter->second, std::map()) );
name = Particle::findFreeName(periode, elementno);
_atom_particlenames[graphiter->second][elementno] = name;
LOG(1, "INFO: Adding particle " << name << " for atom "
<< *walker << " with element " << elementno << " and initial charge "
<< charge << ".");
factory.createInstance(name, elementno, charge);
#ifndef NDEBUG
std::pair< detail::average_charge_map_t::iterator, bool> inserter =
#endif
average_charge_map.insert( std::make_pair(name, std::make_pair(charge,1)) );
ASSERT( inserter.second,
"addToParticleRegistry() - name "+name+" already present in average charge map?");
}
}
// go through map, average and insert average charge into factory
for (detail::average_charge_map_t::const_iterator averageiter = average_charge_map.begin();
averageiter != average_charge_map.end(); ++averageiter) {
const std::string &name = averageiter->first;
const double &charge = averageiter->second.first;
const size_t × = averageiter->second.second;
ASSERT( times != 0,
"addToParticleRegistry() - zero times in average charge map is not allowed.");
Particle *particle = registry.getByName(name);
particle->charge = charge/(double)times;
}
}
static bool isNotHydrogen(const atom * const _atom)
{
return (_atom->getElementNo() != (atomicNumber_t) 1);
}
static struct KeySetSizeComp {
bool operator() (const KeySet &a, const KeySet &b) { return a.size()
#include
#include "Fragmentation/Summation/AllLevelOrthogonalSummator.hpp"
#include "Fragmentation/Summation/IndexSet.hpp"
#include "Fragmentation/Summation/IndexSetContainer.hpp"
#include "Fragmentation/Summation/SubsetMap.hpp"
#include "Fragmentation/Summation/Containers/PartialChargesFused.hpp"
#include "Fragmentation/Summation/Containers/PartialChargesMap.hpp"
#include "Fragmentation/Summation/SetValues/IndexedPartialCharges.hpp"
ActionState::ptr PotentialFitFragmentPartialChargesAction::performCall()
{
// check for selected atoms
const World &world = World::getConstInstance();
const std::vector selected_atoms = world.getSelectedAtoms();
if (selected_atoms.empty()) {
STATUS("There are no atoms selected for fitting partial charges to.");
return Action::failure;
}
LOG(3, "There are " << selected_atoms.size() << " selected atoms.");
/// obtain possible fragments to each selected atom
const AtomFragmentsMap &atomfragments = AtomFragmentsMap::getConstInstance();
if (!atomfragments.checkCompleteness()) {
ELOG(0, "AtomFragmentsMap failed internal consistency check, missing forcekeysets?");
return Action::failure;
}
const std::set fragments =
accumulateKeySetsForAtoms( atomfragments.getMap(), selected_atoms);
const size_t NoNonHydrogens =
std::count_if(selected_atoms.begin(), selected_atoms.end(), isNotHydrogen);
if (fragments.size() < NoNonHydrogens) {
ELOG(0, "Obtained fewer fragments than there are atoms, has AtomFragments been loaded?");
return Action::failure;
}
// get results from container
FragmentationResultContainer& resultscontainer = FragmentationResultContainer::getInstance();
#if defined(HAVE_JOBMARKET) && defined(HAVE_VMG)
const bool DoLongrange = resultscontainer.areFullRangeResultsPresent();
if (!DoLongrange) {
ELOG(1, "FragmentationResultContainer does not contain long-range results, cannot fit partial charges.");
return Action::failure;
}
if (resultscontainer.getKeySets().KeySets.empty()) {
STATUS("There are no results in the container.");
return Action::failure;
}
// reduce given fragments to homologous graphs for later lookup
detail::KeysetsToGraph_t keyset_graphs;
detail::fittedcharges_per_fragment_t fittedcharges_per_fragment;
getKeySetsToGraphMapping(keyset_graphs, fragments, atomfragments);
/// then go through all fragments and get partial charges for each
const bool status = getPartialChargesForAllIndexSets(
fittedcharges_per_fragment,
resultscontainer,
fragments,
params.radius.get(),
params.enforceZeroCharge.get());
if (!status)
return Action::failure;
/// obtain average charge for each atom the fitted charges over all its fragments
detail::fitted_charges_t fitted_charges;
{
// place all into boost::fusion map for summation
std::map PartialCharges_fused;
std::vector indexsets;
JobId_t counter = 0;
for (detail::fittedcharges_per_fragment_t::const_iterator iter = fittedcharges_per_fragment.begin();
iter != fittedcharges_per_fragment.end(); ++iter) {
const KeySet ¤tset = iter->first;
// place as IndexSet into container
{
IndexSet tempset;
tempset.insert(currentset.begin(), currentset.end());
indexsets.push_back(tempset);
}
LOG(3, "Inserting " << iter->first << " with charges " << iter->second << " into chargemap instance.");
PartialChargesMap_t chargemap;
const AtomFragmentsMap::indices_t &full_currentset = atomfragments.getFullKeyset(currentset);
IndexedPartialCharges::indices_t indices(full_currentset.begin(), full_currentset.end());
boost::fusion::at_key(chargemap) =
IndexedPartialCharges(full_currentset, iter->second);
PartialCharges_fused.insert( std::make_pair(counter++, chargemap) );
}
ASSERT( counter == fragments.size(),
"PotentialFitPartialChargesAction::performCall() - not all fragments' keysets were created?");
ASSERT( indexsets.size() == fragments.size(),
"PotentialFitPartialChargesAction::performCall() - not all fragments' keysets were created?");
// prepare index set hierarchy
std::map< JobId_t, size_t > IdentityLookup;
size_t MaxLevel = std::max_element(fragments.begin(), fragments.end(), keyset_comparator)->size();
LOG(3, "The maximum level is " << MaxLevel);
size_t max_indices = fragments.size();
LOG(3, "There are " << max_indices << " keysets for the selected atoms.");
{
for (size_t index = 0; index < max_indices; ++index)
IdentityLookup.insert( std::make_pair( (JobId_t)index, index ) );
}
IndexSetContainer::ptr container(new IndexSetContainer(indexsets));
SubsetMap::ptr subsetmap(new SubsetMap(*container));
// and sum up
PartialChargesMap_t ZeroInstances;
ZeroInstanceInitializer initZeroInstance(ZeroInstances);
boost::mpl::for_each(boost::ref(initZeroInstance));
//!> results per level of summed up partial charges
std::vector Result_PartialCharges_fused(MaxLevel);
//!> results per index set in terms of value and contribution
std::map<
IndexSet::ptr,
std::pair > Result_perIndexSet_PartialCharges;
AllLevelOrthogonalSummator partialchargeSummer(
subsetmap,
PartialCharges_fused,
container->getContainer(),
IdentityLookup,
Result_PartialCharges_fused,
Result_perIndexSet_PartialCharges,
ZeroInstances);
boost::mpl::for_each(boost::ref(partialchargeSummer));
// TODO: place results into fitted_charges
const IndexedPartialCharges::indexedvalues_t indexed_partial_charges =
boost::fusion::at_key(
Result_PartialCharges_fused.back()
).getValues();
for (IndexedPartialCharges::indexedvalues_t::const_iterator iter = indexed_partial_charges.begin();
iter != indexed_partial_charges.end(); ++iter)
fitted_charges.insert( std::make_pair( iter->first, iter->second.charge));
LOG(3, "Summation has brought forth the following charges per atom index: " << fitted_charges);
}
// make a unique list of the present HomologyGraphs
typedef std::set unique_graphs_t;
unique_graphs_t unique_graphs;
for (detail::KeysetsToGraph_t::const_iterator keysetiter = keyset_graphs.begin();
keysetiter != keyset_graphs.end(); ++keysetiter)
unique_graphs.insert(keysetiter->second);
/// make Particles be used for every atom that was fitted on the same number of graphs
detail::GraphIndex_t GraphIndex;
size_t index = 0;
for (unique_graphs_t::const_iterator graphiter = unique_graphs.begin();
graphiter != unique_graphs.end(); ++graphiter) {
GraphIndex.insert( std::make_pair( *graphiter, index++));
}
LOG(2, "DEBUG: There are " << index << " unique graphs in the homology container.");
ASSERT( index == unique_graphs.size(),
"PotentialFitFragmentPartialChargesAction::performCall() - could not enumerate all unique graphs?");
// go through every non-hydrogen atom, get all graphs, convert to GraphIndex and store
detail::GraphIndices_t GraphIndices;
const AtomFragmentsMap::AtomFragmentsMap_t &atommap = atomfragments.getMap();
for (std::vector::const_iterator atomiter = selected_atoms.begin();
atomiter != selected_atoms.end(); ++atomiter) {
// use the non-hydrogen here
const atomId_t walkerid = (*atomiter)->getId();
const atomId_t surrogateid = getNonHydrogenSurrogate(*atomiter)->getId();
if (surrogateid != walkerid)
continue;
const AtomFragmentsMap::AtomFragmentsMap_t::const_iterator keysetsiter =
atommap.find(walkerid);
ASSERT(keysetsiter != atommap.end(),
"PotentialFitFragmentPartialChargesAction::performCall() - we checked already that "
+toString(surrogateid)+" should be present!");
const AtomFragmentsMap::keysets_t & keysets = keysetsiter->second;
// go over all fragments associated to this atom
detail::AtomsGraphIndices_t AtomsGraphIndices;
for (AtomFragmentsMap::keysets_t::const_iterator keysetsiter = keysets.begin();
keysetsiter != keysets.end(); ++keysetsiter) {
const KeySet &keyset = *keysetsiter;
const std::map::const_iterator keysetgraphiter =
keyset_graphs.find(keyset);
ASSERT( keysetgraphiter != keyset_graphs.end(),
"PotentialFitFragmentPartialChargesAction::performCall() - keyset "+toString(keyset)
+" not contained in keyset_graphs.");
const HomologyGraph &graph = keysetgraphiter->second;
const detail::GraphIndex_t::const_iterator indexiter = GraphIndex.find(graph);
ASSERT( indexiter != GraphIndex.end(),
"PotentialFitFragmentPartialChargesAction::performCall() - graph "+toString(graph)
+" not contained in GraphIndex.");
AtomsGraphIndices.insert( indexiter->second );
}
GraphIndices.insert( detail::GraphIndices_t::value_type(AtomsGraphIndices, walkerid) );
LOG(2, "DEBUG: Atom #" << walkerid << "," << **atomiter << ". has graph indices "
<< AtomsGraphIndices);
}
// then graphs from non-hydrogen bond partner for all hydrogens
for (std::vector::const_iterator atomiter = selected_atoms.begin();
atomiter != selected_atoms.end(); ++atomiter) {
// use the non-hydrogen here
const atomId_t walkerid = (*atomiter)->getId();
const atomId_t surrogateid = getNonHydrogenSurrogate((*atomiter))->getId();
if (surrogateid == walkerid)
continue;
detail::GraphIndices_t::right_const_iterator graphiter = GraphIndices.right.find(surrogateid);
ASSERT( graphiter != GraphIndices.right.end(),
"PotentialFitFragmentPartialChargesAction::performCall() - atom #"+toString(surrogateid)
+" not contained in GraphIndices.");
const detail::AtomsGraphIndices_t &AtomsGraphIndices = graphiter->second;
GraphIndices.insert( detail::GraphIndices_t::value_type(AtomsGraphIndices, walkerid) );
LOG(2, "DEBUG: Hydrogen #" << walkerid << ", " << **atomiter
<< ", has graph indices " << AtomsGraphIndices);
}
/// place all fitted (and now averaged) charges into ParticleRegistry
detail::AtomParticleNames_t atom_particlenames;
addToParticleRegistry(
ParticleFactory::getConstInstance(),
ParticleRegistry::getInstance(),
*World::getInstance().getPeriode(),
fitted_charges,
GraphIndices,
atom_particlenames);
for (World::AtomSelectionIterator atomiter = World::getInstance().beginAtomSelection();
atomiter != World::getInstance().endAtomSelection(); ++atomiter) {
atom * const walker = atomiter->second;
const atomId_t walkerid = atomiter->first;
const detail::GraphIndices_t::right_const_iterator graphiter =
GraphIndices.right.find(walkerid);
ASSERT( graphiter != GraphIndices.right.end(),
"PotentialFitFragmentPartialChargesAction::performCall() - cannot find "
+toString(walkerid)+" in GraphIndices.");
const detail::AtomsGraphIndices_t &graphindex = graphiter->second;
const detail::AtomParticleNames_t::const_iterator particlesetiter =
atom_particlenames.find(graphindex);
ASSERT( particlesetiter != atom_particlenames.end(),
"PotentialFitFragmentPartialChargesAction::performCall() - cannot find "
+toString(graphindex)+" in atom_particlenames.");
const std::map::const_iterator nameiter =
particlesetiter->second.find(walker->getElementNo());
ASSERT( nameiter != particlesetiter->second.end(),
"PotentialFitFragmentPartialChargesAction::performCall() - ");
walker->setParticleName(nameiter->second);
LOG(1, "INFO: atom " << *walker << " received the following particle "
<< walker->getParticleName());
}
return Action::success;
#else
ELOG(0, "Long-range support not compiled in, cannot fit partial charges to fragment's potential grids.");
return Action::failure;
#endif
}
ActionState::ptr PotentialFitFragmentPartialChargesAction::performUndo(ActionState::ptr _state) {
return Action::success;
}
ActionState::ptr PotentialFitFragmentPartialChargesAction::performRedo(ActionState::ptr _state){
return Action::success;
}
bool PotentialFitFragmentPartialChargesAction::canUndo() {
return false;
}
bool PotentialFitFragmentPartialChargesAction::shouldUndo() {
return false;
}
/** =========== end of function ====================== */