Changeset b8b75d


Ignore:
Timestamp:
Oct 17, 2009, 7:47:58 PM (16 years ago)
Author:
Frederik Heber <heber@…>
Branches:
Action_Thermostats, Add_AtomRandomPerturbation, Add_FitFragmentPartialChargesAction, Add_RotateAroundBondAction, Add_SelectAtomByNameAction, Added_ParseSaveFragmentResults, AddingActions_SaveParseParticleParameters, 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, Candidate_v1.7.0, 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_BoundInBox_CenterInBox_MoleculeActions, 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, ImplicitCharges, IndependentFragmentGrids, IndependentFragmentGrids_IndividualZeroInstances, IndependentFragmentGrids_IntegrationTest, IndependentFragmentGrids_Sole_NN_Calculation, JobMarket_RobustOnKillsSegFaults, JobMarket_StableWorkerPool, JobMarket_unresolvable_hostname_fix, MoreRobust_FragmentAutomation, 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
Children:
266237
Parents:
07051c
Message:

Refactored CreateAdjacencyList().

Location:
src
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • src/molecule.hpp

    r07051c rb8b75d  
    248248  void CreateAdjacencyList(ofstream *out, double bonddistance, bool IsAngstroem);
    249249  void CreateListOfBondsPerAtom(ofstream *out);
     250  int CorrectBondDegree(ofstream *out);
     251  void OutputBondsList(ofstream *out);
     252  int CountAtomsBonds(int nr);
     253
    250254
    251255  // Graph analysis
  • src/molecule_graph.cpp

    r07051c rb8b75d  
    1111#include "element.hpp"
    1212#include "helpers.hpp"
     13#include "linkedcell.hpp"
    1314#include "lists.hpp"
    1415#include "memoryallocator.hpp"
     
    7374{
    7475
    75   atom *Walker = NULL, *OtherWalker = NULL, *Candidate = NULL;
    76   int No, NoBonds, CandidateBondNo;
    77   int NumberCells, divisor[NDIM], n[NDIM], N[NDIM], index, Index, j;
    78   molecule **CellList;
     76  atom *Walker = NULL;
     77  atom *OtherWalker = NULL;
     78  atom **AtomMap = NULL;
     79  int n[NDIM];
    7980  double distance, MinDistance, MaxDistance;
    80   double *matrix = ReturnFullMatrixforSymmetric(cell_size);
    81   Vector x;
    82   int FalseBondDegree = 0;
     81  LinkedCell *LC = NULL;
     82  LinkedNodes *List = NULL;
     83  LinkedNodes *OtherList = NULL;
    8384
    8485  BondDistance = bonddistance; // * ((IsAngstroem) ? 1. : 1./AtomicLengthToAngstroem);
     
    9495
    9596  if (AtomCount != 0) {
    96     // 1. find divisor for each axis, such that a sphere with radius of at least bonddistance can be placed into each cell
    97     j=-1;
    98     for (int i=0;i<NDIM;i++) {
    99       j += i+1;
    100       divisor[i] = (int)floor(cell_size[j]/bonddistance); // take smaller value such that size of linked cell is at least bonddistance
    101       //*out << Verbose(1) << "divisor[" << i << "]  = " << divisor[i] << "." << endl;
    102     }
    103     // 2a. allocate memory for the cell list
    104     NumberCells = divisor[0]*divisor[1]*divisor[2];
    105     *out << Verbose(1) << "Allocating " << NumberCells << " cells." << endl;
    106     CellList = Malloc<molecule*>(NumberCells, "molecule::CreateAdjacencyList - ** CellList");
    107     for (int i=NumberCells;i--;)
    108       CellList[i] = NULL;
    109 
    110     // 2b. put all atoms into its corresponding list
     97    LC = new LinkedCell(this, bonddistance);
     98
     99    // create a list to map Tesselpoint::nr to atom *
     100    AtomMap = Malloc<atom *>(AtomCount, "molecule::CreateAdjacencyList - **AtomCount");
    111101    Walker = start;
    112     while(Walker->next != end) {
     102    while (Walker->next != end) {
    113103      Walker = Walker->next;
    114       //*out << Verbose(1) << "Current atom is " << *Walker << " with coordinates ";
    115       //Walker->x.Output(out);
    116       //*out << "." << endl;
    117       // compute the cell by the atom's coordinates
    118       j=-1;
    119       for (int i=0;i<NDIM;i++) {
    120         j += i+1;
    121         x.CopyVector(&(Walker->x));
    122         x.KeepPeriodic(out, matrix);
    123         n[i] = (int)floor(x.x[i]/cell_size[j]*(double)divisor[i]);
    124       }
    125       index = n[2] + (n[1] + n[0] * divisor[1]) * divisor[2];
    126       //*out << Verbose(1) << "Atom " << *Walker << " goes into cell number [" << n[0] << "," << n[1] << "," << n[2] << "] = " << index << "." << endl;
    127       // add copy atom to this cell
    128       if (CellList[index] == NULL)  // allocate molecule if not done
    129         CellList[index] = new molecule(elemente);
    130       OtherWalker = CellList[index]->AddCopyAtom(Walker); // add a copy of walker to this atom, father will be walker for later reference
    131       //*out << Verbose(1) << "Copy Atom is " << *OtherWalker << "." << endl;
    132     }
    133     //for (int i=0;i<NumberCells;i++)
    134       //*out << Verbose(1) << "Cell number " << i << ": " << CellList[i] << "." << endl;
    135 
     104      AtomMap[Walker->nr] = Walker;
     105    }
    136106
    137107    // 3a. go through every cell
    138     for (N[0]=divisor[0];N[0]--;)
    139       for (N[1]=divisor[1];N[1]--;)
    140         for (N[2]=divisor[2];N[2]--;) {
    141           Index = N[2] + (N[1] + N[0] * divisor[1]) * divisor[2];
    142           if (CellList[Index] != NULL) { // if there atoms in this cell
    143             //*out << Verbose(1) << "Current cell is " << Index << "." << endl;
    144             // 3b. for every atom therein
    145             Walker = CellList[Index]->start;
    146             while (Walker->next != CellList[Index]->end) {  // go through every atom
    147               Walker = Walker->next;
     108    for (LC->n[0] = 0; LC->n[0] < LC->N[0]; LC->n[0]++)
     109      for (LC->n[1] = 0; LC->n[1] < LC->N[1]; LC->n[1]++)
     110        for (LC->n[2] = 0; LC->n[2] < LC->N[2]; LC->n[2]++) {
     111          List = LC->GetCurrentCell();
     112          //*out << Verbose(2) << "Current cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << " containing " << List->size() << " points." << endl;
     113          if (List != NULL) {
     114            for (LinkedNodes::iterator Runner = List->begin(); Runner != List->end(); Runner++) {
     115              Walker = AtomMap[(*Runner)->nr];
    148116              //*out << Verbose(0) << "Current Atom is " << *Walker << "." << endl;
    149117              // 3c. check for possible bond between each atom in this and every one in the 27 cells
     
    151119                for (n[1]=-1;n[1]<=1;n[1]++)
    152120                  for (n[2]=-1;n[2]<=1;n[2]++) {
    153                      // compute the index of this comparison cell and make it periodic
    154                     index = ((N[2]+n[2]+divisor[2])%divisor[2]) + (((N[1]+n[1]+divisor[1])%divisor[1]) + ((N[0]+n[0]+divisor[0])%divisor[0]) * divisor[1]) * divisor[2];
    155                     //*out << Verbose(1) << "Number of comparison cell is " << index << "." << endl;
    156                     if (CellList[index] != NULL) {  // if there are any atoms in this cell
    157                       OtherWalker = CellList[index]->start;
    158                       while(OtherWalker->next != CellList[index]->end) {  // go through every atom in this cell
    159                         OtherWalker = OtherWalker->next;
    160                         //*out << Verbose(0) << "Current comparison atom is " << *OtherWalker << "." << endl;
    161                         /// \todo periodic check is missing here!
    162                         //*out << Verbose(1) << "Checking distance " << OtherWalker->x.PeriodicDistanceSquared(&(Walker->x), cell_size) << " against typical bond length of " << bonddistance*bonddistance << "." << endl;
    163                         MinDistance = OtherWalker->type->CovalentRadius + Walker->type->CovalentRadius;
    164                         MinDistance *= (IsAngstroem) ? 1. : 1./AtomicLengthToAngstroem;
    165                         MaxDistance = MinDistance + BONDTHRESHOLD;
    166                         MinDistance -= BONDTHRESHOLD;
    167                         distance = OtherWalker->x.PeriodicDistanceSquared(&(Walker->x), cell_size);
    168                         if ((OtherWalker->father->nr > Walker->father->nr) && (distance <= MaxDistance*MaxDistance) && (distance >= MinDistance*MinDistance)) { // create bond if distance is smaller
    169                           //*out << Verbose(1) << "Adding Bond between " << *Walker << " and " << *OtherWalker << " in distance " << sqrt(distance) << "." << endl;
    170                           AddBond(Walker->father, OtherWalker->father, 1);  // also increases molecule::BondCount
    171                         } else {
    172                           //*out << Verbose(1) << "Not Adding: Wrong label order or distance too great." << endl;
     121                    OtherList = LC->GetRelativeToCurrentCell(n);
     122                    //*out << Verbose(2) << "Current relative cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << " containing " << List->size() << " points." << endl;
     123                    if (OtherList != NULL) {
     124                      for (LinkedNodes::iterator OtherRunner = OtherList->begin(); OtherRunner != OtherList->end(); OtherRunner++) {
     125                        if ((*OtherRunner)->nr > Walker->nr) {
     126                          OtherWalker = AtomMap[(*OtherRunner)->nr];
     127                          //*out << Verbose(1) << "Checking distance " << OtherWalker->x.PeriodicDistanceSquared(&(Walker->x), cell_size) << " against typical bond length of " << bonddistance*bonddistance << "." << endl;
     128                          MinDistance = OtherWalker->type->CovalentRadius + Walker->type->CovalentRadius;
     129                          MinDistance *= (IsAngstroem) ? 1. : 1./AtomicLengthToAngstroem;
     130                          MaxDistance = MinDistance + BONDTHRESHOLD;
     131                          MinDistance -= BONDTHRESHOLD;
     132                          distance = OtherWalker->x.PeriodicDistanceSquared(&(Walker->x), cell_size);
     133                          if ((OtherWalker->father->nr > Walker->father->nr) && (distance <= MaxDistance*MaxDistance) && (distance >= MinDistance*MinDistance)) { // create bond if distance is smaller
     134                            //*out << Verbose(1) << "Adding Bond between " << *Walker << " and " << *OtherWalker << " in distance " << sqrt(distance) << "." << endl;
     135                            AddBond(Walker->father, OtherWalker->father, 1);  // also increases molecule::BondCount
     136                          } else {
     137                            //*out << Verbose(1) << "Not Adding: Wrong label order or distance too great." << endl;
     138                          }
    173139                        }
    174140                      }
     
    178144          }
    179145        }
    180 
    181 
    182 
    183     // 4. free the cell again
    184     for (int i=NumberCells;i--;)
    185       if (CellList[i] != NULL) {
    186         delete(CellList[i]);
    187       }
    188     Free(&CellList);
     146    Free(&AtomMap);
     147    delete(LC);
     148    *out << Verbose(1) << "I detected " << BondCount << " bonds in the molecule with distance " << BondDistance << "." << endl;
    189149
    190150    // create the adjacency list per atom
    191151    CreateListOfBondsPerAtom(out);
    192152
    193     // correct Bond degree of each bond by checking both bond partners for a mismatch between valence and current sum of bond degrees,
    194     // iteratively increase the one first where the other bond partner has the fewest number of bonds (i.e. in general bonds oxygene
    195     // preferred over carbon bonds). Beforehand, we had picked the first mismatching partner, which lead to oxygenes with single instead of
    196     // double bonds as was expected.
    197     if (BondCount != 0) {
    198       NoCyclicBonds = 0;
    199       *out << Verbose(1) << "Correcting Bond degree of each bond ... ";
    200       do {
    201         No = 0; // No acts as breakup flag (if 1 we still continue)
    202         Walker = start;
    203         while (Walker->next != end) { // go through every atom
    204           Walker = Walker->next;
    205           // count valence of first partner
    206           NoBonds = 0;
    207           for(j=0;j<NumberOfBondsPerAtom[Walker->nr];j++)
    208             NoBonds += ListOfBondsPerAtom[Walker->nr][j]->BondDegree;
    209           *out << Verbose(3) << "Walker " << *Walker << ": " << (int)Walker->type->NoValenceOrbitals << " > " << NoBonds << "?" << endl;
    210           if ((int)(Walker->type->NoValenceOrbitals) > NoBonds) { // we have a mismatch, check all bonding partners for mismatch
    211             Candidate = NULL;
    212             CandidateBondNo = -1;
    213             for(int i=0;i<NumberOfBondsPerAtom[Walker->nr];i++) { // go through each of its bond partners
    214               OtherWalker = ListOfBondsPerAtom[Walker->nr][i]->GetOtherAtom(Walker);
    215               // count valence of second partner
    216               NoBonds = 0;
    217               for(j=0;j<NumberOfBondsPerAtom[OtherWalker->nr];j++)
    218                 NoBonds += ListOfBondsPerAtom[OtherWalker->nr][j]->BondDegree;
    219               *out << Verbose(3) << "OtherWalker " << *OtherWalker << ": " << (int)OtherWalker->type->NoValenceOrbitals << " > " << NoBonds << "?" << endl;
    220               if ((int)(OtherWalker->type->NoValenceOrbitals) > NoBonds) { // check if possible candidate
    221                 if ((Candidate == NULL) || (NumberOfBondsPerAtom[Candidate->nr] > NumberOfBondsPerAtom[OtherWalker->nr])) { // pick the one with fewer number of bonds first
    222                   Candidate = OtherWalker;
    223                   CandidateBondNo = i;
    224                   *out << Verbose(3) << "New candidate is " << *Candidate << "." << endl;
    225                 }
    226               }
    227             }
    228             if ((Candidate != NULL) && (CandidateBondNo != -1)) {
    229               ListOfBondsPerAtom[Walker->nr][CandidateBondNo]->BondDegree++;
    230               *out << Verbose(2) << "Increased bond degree for bond " << *ListOfBondsPerAtom[Walker->nr][CandidateBondNo] << "." << endl;
    231             } else
    232               *out << Verbose(2) << "Could not find correct degree for atom " << *Walker << "." << endl;
    233               FalseBondDegree++;
    234           }
    235         }
    236       } while (No);
    237     *out << " done." << endl;
    238     } else
    239       *out << Verbose(1) << "BondCount is " << BondCount << ", no bonds between any of the " << AtomCount << " atoms." << endl;
    240     *out << Verbose(1) << "I detected " << BondCount << " bonds in the molecule with distance " << bonddistance << ", " << FalseBondDegree << " bonds could not be corrected." << endl;
     153    // correct bond degree by comparing valence and bond degree
     154    CorrectBondDegree(out);
    241155
    242156    // output bonds for debugging (if bond chain list was correctly installed)
    243     *out << Verbose(1) << endl << "From contents of bond chain list:";
    244     bond *Binder = first;
    245     while(Binder->next != last) {
    246       Binder = Binder->next;
    247       *out << *Binder << "\t" << endl;
    248     }
    249     *out << endl;
     157    OutputBondsList(out);
    250158  } else
    251159    *out << Verbose(1) << "AtomCount is " << AtomCount << ", thus no bonds, no connections!." << endl;
    252160  *out << Verbose(0) << "End of CreateAdjacencyList." << endl;
    253   Free(&matrix);
    254 
     161};
     162
     163/** Prints a list of all bonds to \a *out.
     164 * \param output stream
     165 */
     166void molecule::OutputBondsList(ofstream *out)
     167{
     168  *out << Verbose(1) << endl << "From contents of bond chain list:";
     169  bond *Binder = first;
     170  while(Binder->next != last) {
     171    Binder = Binder->next;
     172    *out << *Binder << "\t" << endl;
     173  }
     174  *out << endl;
     175};
     176
     177/** Sums up the number of bonds times bond::BondDegree the atom with index \a nr has.
     178 * \param nr indexof atom
     179 * \return number of bonds each weighted with its bond::BondDegree
     180 */
     181int molecule::CountAtomsBonds(int nr)
     182{
     183  int NoBonds = 0;
     184  for(int j=0;j<NumberOfBondsPerAtom[nr];j++)
     185    NoBonds += ListOfBondsPerAtom[nr][j]->BondDegree;
     186  return NoBonds;
     187};
     188
     189/** correct bond degree by comparing valence and bond degree.
     190 * correct Bond degree of each bond by checking both bond partners for a mismatch between valence and current sum of bond degrees,
     191 * iteratively increase the one first where the other bond partner has the fewest number of bonds (i.e. in general bonds oxygene
     192 * preferred over carbon bonds). Beforehand, we had picked the first mismatching partner, which lead to oxygenes with single instead of
     193 * double bonds as was expected.
     194 * \param *out output stream for debugging
     195 * \return number of bonds that could not be corrected
     196 */
     197int molecule::CorrectBondDegree(ofstream *out)
     198{
     199  int No = 0;
     200  int NoBonds = 0;
     201  int CandidateBondNo = 0;
     202  int FalseBondDegree = 0;
     203  atom *Walker = NULL;
     204  atom *Candidate = NULL;
     205  atom *OtherWalker = NULL;
     206
     207  if (BondCount != 0) {
     208    NoCyclicBonds = 0;
     209    *out << Verbose(1) << "Correcting Bond degree of each bond ... ";
     210    do {
     211      No = 0; // No acts as breakup flag (if 1 we still continue)
     212      Walker = start;
     213      while (Walker->next != end) { // go through every atom
     214        Walker = Walker->next;
     215
     216        NoBonds = CountAtomsBonds(Walker->nr);
     217        *out << Verbose(3) << "Walker " << *Walker << ": " << (int)Walker->type->NoValenceOrbitals << " > " << NoBonds << "?" << endl;
     218        if ((int)(Walker->type->NoValenceOrbitals) > NoBonds) { // we have a mismatch, check all bonding partners for mismatch
     219          Candidate = NULL;
     220          CandidateBondNo = -1;
     221          for(int i=0;i<NumberOfBondsPerAtom[Walker->nr];i++) { // go through each of its bond partners
     222            OtherWalker = ListOfBondsPerAtom[Walker->nr][i]->GetOtherAtom(Walker);
     223            NoBonds = CountAtomsBonds(OtherWalker->nr);
     224            *out << Verbose(3) << "OtherWalker " << *OtherWalker << ": " << (int)OtherWalker->type->NoValenceOrbitals << " > " << NoBonds << "?" << endl;
     225            if ((int)(OtherWalker->type->NoValenceOrbitals) > NoBonds) { // check if possible candidate
     226              if ((Candidate == NULL) || (NumberOfBondsPerAtom[Candidate->nr] > NumberOfBondsPerAtom[OtherWalker->nr])) { // pick the one with fewer number of bonds first
     227                Candidate = OtherWalker;
     228                CandidateBondNo = i;
     229                *out << Verbose(3) << "New candidate is " << *Candidate << "." << endl;
     230              }
     231            }
     232          }
     233          if ((Candidate != NULL) && (CandidateBondNo != -1)) {
     234            ListOfBondsPerAtom[Walker->nr][CandidateBondNo]->BondDegree++;
     235            *out << Verbose(2) << "Increased bond degree for bond " << *ListOfBondsPerAtom[Walker->nr][CandidateBondNo] << "." << endl;
     236          } else {
     237            *out << Verbose(2) << "Could not find correct degree for atom " << *Walker << "." << endl;
     238            FalseBondDegree++;
     239          }
     240        }
     241      }
     242    } while (No);
     243    *out << " done." << endl;
     244  } else {
     245    *out << Verbose(1) << "BondCount is " << BondCount << ", no bonds between any of the " << AtomCount << " atoms." << endl;
     246  }
     247  *out << FalseBondDegree << " bonds could not be corrected." << endl;
     248
     249  return (FalseBondDegree);
    255250};
    256251
     
    285280  return No;
    286281};
     282
     283
    287284/** Returns Shading as a char string.
    288285 * \param color the Shading
Note: See TracChangeset for help on using the changeset viewer.