- fixed issues that popped up after clean-up, mostly namespaces
authorAnatoly <Titov_AI@pnpi.nrcki.ru>
Wed, 21 Oct 2020 10:20:54 +0000 (13:20 +0300)
committerAnatoly <Titov_AI@pnpi.nrcki.ru>
Wed, 21 Oct 2020 10:20:54 +0000 (13:20 +0300)
src/colorvec.cpp

index b562ba91a900178192202cab8b6a345aa2d4c944..0b5fd61e84415f385273870224497b377880de26 100644 (file)
@@ -40,9 +40,6 @@
  * \ingroup module_trajectoryanalysis
  */
 
-#include "domaintype.h"
-#include "newfit.h"
-
 #include <gromacs/trajectoryanalysis.h>
 #include <gromacs/trajectoryanalysis/topologyinformation.h>
 #include <gromacs/selection/nbsearch.h>
@@ -54,6 +51,7 @@
 #include <omp.h>
 #include <vector>
 #include <cstdio>
+#include <iomanip>
 
 // структура углов для одной краски
 struct colorLocalAngles {
@@ -121,7 +119,7 @@ void betaListDigestion(const std::string &inputFile, std::vector< std::vector< s
 }
 
 // функция выделения индексов для бэталистов
-inline void aminoacidsIndexation(const std::vector< size_t > inputIndex, const TopologyInformation &top, std::vector< std::vector< size_t > > &aminoacidsIndex) {
+inline void aminoacidsIndexation(const std::vector< size_t > inputIndex, const gmx::TopologyInformation &top, std::vector< std::vector< size_t > > &aminoacidsIndex) {
     aminoacidsIndex.resize(0);
     for (size_t i = 0; i < inputIndex.size(); ++i) {
         aminoacidsIndex.resize(std::max(aminoacidsIndex.size(), static_cast< size_t >(top.atoms()->atom[inputIndex[i]].resind)));
@@ -130,7 +128,7 @@ inline void aminoacidsIndexation(const std::vector< size_t > inputIndex, const T
 }
 
 // функция поиска RVec в кадре по имени->индексу
-gmx::RVec returnRVec(const std::vector< RVec > &frame, const std::vector< std::pair< std::string, size_t > > &colorsIndex, std::string toFind) {
+gmx::RVec returnRVec(const std::vector< gmx::RVec > &frame, const std::vector< std::pair< std::string, size_t > > &colorsIndex, std::string toFind) {
     for (auto &i : colorsIndex) {
         if (i.first == toFind) {
             return frame[i.second];
@@ -152,7 +150,7 @@ inline double RVecAngle(const gmx::RVec &a, const gmx::RVec &b) {
 }
 
 // вычисление внутренних углов в краске
-void colorsAnglesEvaluation(const std::vector< RVec > &frame, const std::vector< std::vector< std::pair< std::string, size_t > > > &colorsIndex, std::vector< colorLocalAngles > &colorStruct) {
+void colorsAnglesEvaluation(const std::vector< gmx::RVec > &frame, const std::vector< std::vector< std::pair< std::string, size_t > > > &colorsIndex, std::vector< colorLocalAngles > &colorStruct) {
     colorStruct.resize(colorsIndex.size());
     #pragma omp parallel for ordered schedule(dynamic)
     for (size_t i = 0; i < colorsIndex.size(); ++i) {
@@ -180,7 +178,7 @@ void colorsAnglesEvaluation(const std::vector< RVec > &frame, const std::vector<
 }
 
 // вычисление направляющего вектора в бэта-листе
-inline void betaListsRVecsEvaluation(const std::vector< RVec > &frame, const std::vector< std::vector< unsigned int > > &inputBetaLists, std::vector< gmx::RVec > &temp, const std::vector< size_t > &inputCA) {
+inline void betaListsRVecsEvaluation(const std::vector< gmx::RVec > &frame, const std::vector< std::vector< unsigned int > > &inputBetaLists, std::vector< gmx::RVec > &temp, const std::vector< size_t > &inputCA) {
     temp.resize(0);
     gmx::RVec tempA;
     for (const auto &i : inputBetaLists) {
@@ -191,7 +189,7 @@ inline void betaListsRVecsEvaluation(const std::vector< RVec > &frame, const std
 }
 
 // определение близко ли к белку находится краска
-bool isNearPeptide(const t_trxframe &fr, const t_pbc *pbc, const std::vector< RVec > inputFrame,
+bool isNearPeptide(const t_trxframe &fr, const t_pbc *pbc, const std::vector< gmx::RVec > inputFrame,
                           gmx::AnalysisNeighborhood &nbhood, const std::vector< size_t > &inputIndex,
                           const std::vector< std::pair< std::string, size_t > > &inputColor) {
     gmx::AnalysisNeighborhoodSearch      nbsearch = nbhood.initSearch(pbc, gmx::AnalysisNeighborhoodPositions(fr.x, fr.natoms));
@@ -209,7 +207,7 @@ bool isNearPeptide(const t_trxframe &fr, const t_pbc *pbc, const std::vector< RV
 }
 
 // поиск ближайших к краске бэта-листов
-inline void searchNearBetaLists(const t_trxframe &fr, const t_pbc *pbc, const std::vector< RVec > inputFrame,
+inline void searchNearBetaLists(const t_trxframe &fr, const t_pbc *pbc, const std::vector< gmx::RVec > inputFrame,
                                 gmx::AnalysisNeighborhood &nbhood,
                                 const std::vector< std::vector< unsigned int > > inputBLists,
                                 const std::vector< std::pair< std::string, size_t > > &inputColor,
@@ -285,7 +283,7 @@ void anglesFileDump(const int frameNum, const std::string &output, const std::ve
 /*! \brief
  * \ingroup module_trajectoryanalysis
  */
-class colorVec : public TrajectoryAnalysisModule
+class colorVec : public gmx::TrajectoryAnalysisModule
 {
     public:
 
@@ -293,18 +291,18 @@ class colorVec : public TrajectoryAnalysisModule
         virtual ~colorVec();
 
         //! Set the options and setting
-        virtual void initOptions(IOptionsContainer          *options,
-                                 TrajectoryAnalysisSettings *settings);
+        virtual void initOptions(gmx::IOptionsContainer          *options,
+                                 gmx::TrajectoryAnalysisSettings *settings);
 
         //! First routine called by the analysis framework
         // virtual void initAnalysis(const t_trxframe &fr, t_pbc *pbc);
-        virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
-                                  const TopologyInformation        &top);
+        virtual void initAnalysis(const gmx::TrajectoryAnalysisSettings &settings,
+                                  const gmx::TopologyInformation        &top);
 
         //! Call for each frame of the trajectory
         // virtual void analyzeFrame(const t_trxframe &fr, t_pbc *pbc);
         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
-                                  TrajectoryAnalysisModuleData *pdata);
+                                  gmx::TrajectoryAnalysisModuleData *pdata);
 
         //! Last routine called by the analysis framework
         // virtual void finishAnalysis(t_pbc *pbc);
@@ -315,7 +313,7 @@ class colorVec : public TrajectoryAnalysisModule
 
     private:
 
-        SelectionList                                                   sel_;
+        gmx::SelectionList                                              sel_;
         std::string                                                     fnOut;                  // selectable
         std::string                                                     fnBetaListsDat;         // selectable
         float                                                           effRad          {0.8};  // selectable
@@ -324,7 +322,7 @@ class colorVec : public TrajectoryAnalysisModule
         std::vector< std::vector< size_t > >                            aminoacidsIndex;
         std::vector< std::vector< std::pair< std::string, size_t > > >  colorsNames;
         std::vector< std::vector< std::vector< unsigned int > > >       betaLists;
-        AnalysisNeighborhood                                            nb_;
+        gmx::AnalysisNeighborhood                                       nb_;
 
         // Copy and assign disallowed by base.
 };
@@ -346,8 +344,8 @@ colorVec::~colorVec()
  */
 
 void
-colorVec::initOptions(IOptionsContainer         *options,
-                     TrajectoryAnalysisSettings *settings)
+colorVec::initOptions(  gmx::IOptionsContainer          *options,
+                        gmx::TrajectoryAnalysisSettings *settings)
 {
     static const char *const desc[] = {
         "[THISMODULE] to be done"
@@ -355,46 +353,46 @@ colorVec::initOptions(IOptionsContainer         *options,
     // Add the descriptive text (program help text) to the options
     settings->setHelpText(desc);
     // Add option for selection list
-    options->addOption(SelectionOption("sel")
+    options->addOption(gmx::SelectionOption("sel")
                             .storeVector(&sel_)
                             .required().dynamicMask().multiValue()
                             .description("select residue and domains to form a rigid skeleton / -sf"));
     // Add option for output file name
-    options->addOption(FileNameOption("dat").filetype(eftOptionFileType_NR).inputFile()
+    options->addOption(gmx::FileNameOption("dat").filetype(gmx::eftOptionFileType_NR).inputFile()
                             .store(&fnBetaListsDat).description("a file to make dynamic beta lists"));
     // Add option for output file name
-    options->addOption(FileNameOption("out").filetype(eftOptionFileType_NR).outputFile()
+    options->addOption(gmx::FileNameOption("out").filetype(gmx::eftOptionFileType_NR).outputFile()
                             .store(&fnOut).defaultBasename("colorAngles")
                             .description("Index file for the algorithm output."));
     // Add option for effRad constant
-    options->addOption(FloatOption("efRad")
+    options->addOption(gmx::FloatOption("efRad")
                             .store(&effRad)
                             .description("max distance from colors to peptide in nm to consider to be \"near\""));
     // Control input settings
-    settings->setFlags(TrajectoryAnalysisSettings::efNoUserPBC);
-    settings->setFlag(TrajectoryAnalysisSettings::efUseTopX);
+    settings->setFlags(gmx::TrajectoryAnalysisSettings::efNoUserPBC);
+    settings->setFlag(gmx::TrajectoryAnalysisSettings::efUseTopX);
     //settings->setFlag(TrajectoryAnalysisSettings::efRequireTop);
     settings->setPBC(true);
 }
 
 void
-colorVec::initAnalysis(const TrajectoryAnalysisSettings &settings,
-                      const TopologyInformation         &top)
+colorVec::initAnalysis( const gmx::TrajectoryAnalysisSettings   &settings,
+                        const gmx::TopologyInformation          &top)
 {
     // считывание индекса
     index.resize(0);
-    for (ArrayRef< const int >::iterator ai {sel_.front().atomIndices().begin()}; (ai < sel_.front().atomIndices().end()); ai++) {
+    for (gmx::ArrayRef< const int >::iterator ai {sel_.front().atomIndices().begin()}; (ai < sel_.front().atomIndices().end()); ai++) {
         index.push_back(static_cast< size_t >(*ai));
     }
     // считывание индекса
     indexCA.resize(0);
-    for (ArrayRef< const int >::iterator ai {sel_[1].atomIndices().begin()}; (ai < sel_[2].atomIndices().end()); ai++) {
+    for (gmx::ArrayRef< const int >::iterator ai {sel_[1].atomIndices().begin()}; (ai < sel_[2].atomIndices().end()); ai++) {
         index.push_back(static_cast< size_t >(*ai));
     }
     // считывание красок
     colorsNames.resize(sel_.size() - 2);
     for (size_t i = 2; i < sel_.size(); ++i) {
-        for (ArrayRef< const int >::iterator ai {sel_[i].atomIndices().begin()}; (ai < sel_[i].atomIndices().end()); ai++) {
+        for (gmx::ArrayRef< const int >::iterator ai {sel_[i].atomIndices().begin()}; (ai < sel_[i].atomIndices().end()); ai++) {
             colorsNames[i - 2].push_back(std::make_pair(*(top.atoms()->resinfo[top.atoms()->atom[*ai].resind].name), static_cast< size_t >(*ai)));
         }
     }
@@ -416,9 +414,9 @@ colorVec::initAnalysis(const TrajectoryAnalysisSettings &settings,
 }
 
 void
-colorVec::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc, TrajectoryAnalysisModuleData *pdata)
+colorVec::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc, gmx::TrajectoryAnalysisModuleData *pdata)
 {
-    std::vector< RVec > trajectoryFrame;
+    std::vector< gmx::RVec > trajectoryFrame;
     trajectoryFrame.resize(0);
     // считывания текущего фрейма траектории
     for (size_t i {0}; i < index.size(); ++i) {