Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / nbnxm / pairsearch.cpp
blob021f6cb25a4f1044df4167de6198b058a0ddc267
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS 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 GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 /*! \internal \file
37 * \brief
38 * Implements the PairSearch class
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_nbnxm
44 #include "gmxpre.h"
46 #include "pairsearch.h"
48 #include "gromacs/nbnxm/pairlist.h"
49 #include "gromacs/utility/smalloc.h"
52 void SearchCycleCounting::printCycles(FILE *fp,
53 gmx::ArrayRef<const PairsearchWork> work) const
55 fprintf(fp, "\n");
56 fprintf(fp, "ns %4d grid %4.1f search %4.1f",
57 cc_[enbsCCgrid].count(),
58 cc_[enbsCCgrid].averageMCycles(),
59 cc_[enbsCCsearch].averageMCycles());
61 if (work.size() > 1)
63 if (cc_[enbsCCcombine].count() > 0)
65 fprintf(fp, " comb %5.2f",
66 cc_[enbsCCcombine].averageMCycles());
68 fprintf(fp, " s. th");
69 for (const PairsearchWork &workEntry : work)
71 fprintf(fp, " %4.1f",
72 workEntry.cycleCounter.averageMCycles());
75 fprintf(fp, "\n");
78 /*! \brief Frees the contents of a legacy t_nblist struct */
79 static void free_nblist(t_nblist *nl)
81 sfree(nl->iinr);
82 sfree(nl->gid);
83 sfree(nl->shift);
84 sfree(nl->jindex);
85 sfree(nl->jjnr);
86 sfree(nl->excl_fep);
89 #ifndef DOXYGEN
91 PairsearchWork::PairsearchWork() :
92 cp0({{0}}
94 buffer_flags({0, nullptr, 0}),
95 ndistc(0),
96 nbl_fep(new t_nblist),
97 cp1({{0}})
99 nbnxn_init_pairlist_fep(nbl_fep.get());
102 #endif // !DOXYGEN
104 PairsearchWork::~PairsearchWork()
106 sfree(buffer_flags.flag);
108 free_nblist(nbl_fep.get());
111 PairSearch::PairSearch(const int ePBC,
112 const ivec *numDDCells,
113 const gmx_domdec_zones_t *ddZones,
114 const PairlistType pairlistType,
115 const bool haveFep,
116 const int maxNumThreads,
117 gmx::PinningPolicy pinningPolicy) :
118 gridSet_(ePBC, numDDCells, ddZones, pairlistType, haveFep, maxNumThreads, pinningPolicy),
119 work_(maxNumThreads)
121 cycleCounting_.recordCycles_ = (getenv("GMX_NBNXN_CYCLE") != nullptr);