Remove logging from hardware detection
[alexxy/gromacs.git] / api / gmxapi / cpp / context.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2020, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \file
36  * \brief Implementation details of gmxapi::Context
37  *
38  * \todo Share mdrun input handling implementation via modernized modular options framework.
39  * Initial implementation of `launch` relies on borrowed code from the mdrun command
40  * line input processing.
41  *
42  * \author M. Eric Irrgang <ericirrgang@gmail.com>
43  * \ingroup gmxapi
44  */
45
46 #include "gmxapi/context.h"
47
48 #include <cstring>
49
50 #include <memory>
51 #include <utility>
52 #include <vector>
53
54 #include "gromacs/commandline/pargs.h"
55 #include "gromacs/commandline/filenm.h"
56 #include "gromacs/commandline/pargs.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/hardware/detecthardware.h"
59 #include "gromacs/hardware/hw_info.h"
60 #include "gromacs/mdlib/stophandler.h"
61 #include "gromacs/mdrunutility/logging.h"
62 #include "gromacs/mdrunutility/multisim.h"
63 #include "gromacs/mdrun/runner.h"
64 #include "gromacs/mdrunutility/handlerestart.h"
65 #include "gromacs/utility/arraysize.h"
66 #include "gromacs/utility/basenetwork.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/gmxmpi.h"
69 #include "gromacs/utility/init.h"
70 #include "gromacs/utility/physicalnodecommunicator.h"
71
72 #include "gmxapi/mpi/resourceassignment.h"
73 #include "gmxapi/exceptions.h"
74 #include "gmxapi/session.h"
75 #include "gmxapi/status.h"
76 #include "gmxapi/version.h"
77
78 #include "context_impl.h"
79 #include "createsession.h"
80 #include "session_impl.h"
81 #include "workflow.h"
82
83 namespace gmxapi
84 {
85
86 // Support some tag dispatch. Warning: These are just aliases (not strong types).
87 /*!
88  * \brief Logical helper alias to convert preprocessor constant to type.
89  *
90  * \tparam Value Provide the GMX\_LIB\_MPI macro.
91  */
92 template<bool Value>
93 using hasLibraryMpi = std::bool_constant<Value>;
94 /* Note that a no-MPI build still uses the tMPI headers to define MPI_Comm for the
95  * gmx::SimulationContext definition. The dispatching in this file accounts for
96  * these two definitions of SimulationContext. gmxThreadMpi here does not imply
97  * that the library was necessarily compiled with thread-MPI enabled.
98  */
99 using gmxThreadMpi = hasLibraryMpi<false>;
100 using gmxLibMpi    = hasLibraryMpi<true>;
101 using MpiType      = std::conditional_t<GMX_LIB_MPI, gmxLibMpi, gmxThreadMpi>;
102
103 using MpiContextInitializationError = BasicException<struct MpiContextInitialization>;
104
105
106 /*!
107  * \brief Helpers to evaluate the correct precondition for the library build.
108  *
109  * TODO: (#3650) Consider distinct MpiContextManager types for clearer definition of preconditions.
110  */
111 namespace
112 {
113
114 [[maybe_unused]] MPI_Comm validCommunicator(const MPI_Comm& communicator, const gmxThreadMpi&)
115 {
116     if (communicator != MPI_COMM_NULL)
117     {
118         throw MpiContextInitializationError(
119                 "Provided communicator must be MPI_COMM_NULL for GROMACS built without MPI "
120                 "library.");
121     }
122     return communicator;
123 }
124
125 [[maybe_unused]] MPI_Comm validCommunicator(const MPI_Comm& communicator, const gmxLibMpi&)
126 {
127     if (communicator == MPI_COMM_NULL)
128     {
129         throw MpiContextInitializationError("MPI-enabled GROMACS requires a valid communicator.");
130     }
131     return communicator;
132 }
133
134 /*!
135  * \brief Return the communicator if it is appropriate for the environment.
136  *
137  * \throws MpiContextInitializationError if communicator does not match the
138  *  MpiContextManager precondition for the current library configuration.
139  */
140 MPI_Comm validCommunicator(const MPI_Comm& communicator)
141 {
142     return validCommunicator(communicator, MpiType());
143 }
144
145 //! \brief Provide a reasonable default value.
146 MPI_Comm validCommunicator()
147 {
148     return GMX_LIB_MPI ? MPI_COMM_WORLD : MPI_COMM_NULL;
149 }
150
151 } // anonymous namespace
152
153 MpiContextManager::MpiContextManager(MPI_Comm communicator) :
154     communicator_(std::make_unique<MPI_Comm>(validCommunicator(communicator)))
155 {
156     // Safely increments the GROMACS MPI initialization counter after checking
157     // whether the MPI library is already initialized. After this call, MPI_Init
158     // or MPI_Init_thread has been called exactly once.
159     gmx::init(nullptr, nullptr);
160     GMX_RELEASE_ASSERT(!GMX_LIB_MPI || gmx_mpi_initialized(),
161                        "MPI should be initialized before reaching this point.");
162     if (this->communicator() != MPI_COMM_NULL)
163     {
164         // Synchronise at the point of acquiring a MpiContextManager.
165         gmx_barrier(this->communicator());
166     }
167 };
168
169 MpiContextManager::~MpiContextManager()
170 {
171     if (communicator_)
172     {
173         // This is always safe to call. It is a no-op if
174         // thread-MPI, and if the constructor completed then the
175         // MPI library is initialized with reference counting.
176         gmx::finalize();
177     }
178 }
179
180 MpiContextManager::MpiContextManager() : MpiContextManager(validCommunicator()) {}
181
182 MPI_Comm MpiContextManager::communicator() const
183 {
184     if (!communicator_)
185     {
186         throw UsageError("Invalid MpiContextManager. Accessed after `move`?");
187     }
188     return *communicator_;
189 }
190
191 ContextImpl::~ContextImpl() = default;
192
193 [[maybe_unused]] static Context createContext(const ResourceAssignment& resources, const gmxLibMpi&)
194 {
195     CommHandle handle;
196     resources.applyCommunicator(&handle);
197     if (handle.communicator == MPI_COMM_NULL)
198     {
199         throw UsageError("MPI-enabled Simulator contexts require a valid communicator.");
200     }
201     auto contextmanager = MpiContextManager(handle.communicator);
202     auto impl           = ContextImpl::create(std::move(contextmanager));
203     GMX_ASSERT(impl, "ContextImpl creation method should not be able to return null.");
204     auto context = Context(impl);
205     return context;
206 }
207
208 [[maybe_unused]] static Context createContext(const ResourceAssignment& resources, const gmxThreadMpi&)
209 {
210     if (resources.size() > 1)
211     {
212         throw UsageError("Only one thread-MPI Simulation per Context is supported.");
213     }
214     // Thread-MPI Context does not yet have a need for user-provided resources.
215     // However, see #3650.
216     return createContext();
217 }
218
219 Context createContext(const ResourceAssignment& resources)
220 {
221     return createContext(resources, hasLibraryMpi<GMX_LIB_MPI>());
222 }
223
224 Context createContext()
225 {
226     MpiContextManager contextmanager;
227     auto              impl = ContextImpl::create(std::move(contextmanager));
228     GMX_ASSERT(impl, "ContextImpl creation method should not be able to return null.");
229     auto context = Context(impl);
230     return context;
231 }
232
233 ContextImpl::ContextImpl(MpiContextManager&& mpi) noexcept(std::is_nothrow_constructible_v<gmx::LegacyMdrunOptions>) :
234     mpi_(std::move(mpi)),
235     hardwareInformation_(gmx_detect_hardware(
236             gmx::PhysicalNodeCommunicator(mpi_.communicator(), gmx_physicalnode_id_hash())))
237 {
238     // Confirm our understanding of the MpiContextManager invariant.
239     GMX_ASSERT(mpi_.communicator() == MPI_COMM_NULL ? !GMX_LIB_MPI : GMX_LIB_MPI,
240                "Precondition violated: inappropriate communicator for the library environment.");
241     // Make sure we didn't change the data members and overlook implementation details.
242     GMX_ASSERT(session_.expired(),
243                "This implementation assumes an expired weak_ptr at initialization.");
244 }
245
246 std::shared_ptr<ContextImpl> ContextImpl::create(MpiContextManager&& mpi)
247 {
248     std::shared_ptr<ContextImpl> impl;
249     impl.reset(new ContextImpl(std::move(mpi)));
250     return impl;
251 }
252
253 std::shared_ptr<Session> ContextImpl::launch(const Workflow& work)
254 {
255     using namespace gmx;
256     // Much of this implementation is not easily testable: we need tools to inspect simulation
257     // results and to modify simulation inputs.
258
259     std::shared_ptr<Session> launchedSession = nullptr;
260
261     // This implementation can only run one workflow at a time.
262     // Check whether we are already aware of an active session.
263     if (session_.expired())
264     {
265         // Check workflow spec, build graph for current context, launch and return new session.
266         // \todo This is specific to the session implementation...
267         auto        mdNode = work.getNode("MD");
268         std::string filename{};
269         if (mdNode != nullptr)
270         {
271             filename = mdNode->params();
272         }
273
274         /* Mock up the argv interface used by option processing infrastructure.
275          *
276          * As default behavior, automatically extend trajectories from the checkpoint file.
277          * In the future, our API for objects used to initialize a simulation needs to address the fact that currently a
278          * microstate requires data from both the TPR and checkpoint file to be fully specified. Put another way,
279          * current
280          * GROMACS simulations can take a "configuration" as input that does not constitute a complete microstate in
281          * terms of hidden degrees of freedom (integrator/thermostat/barostat/PRNG state), but we want a clear notion of
282          * a microstate for gmxapi interfaces.
283          *
284          * TODO: Remove `-s` and `-cpi` arguments.
285          *       Ref: https://gitlab.com/gromacs/gromacs/-/issues/3652
286          */
287
288         // Set input TPR name
289         mdArgs_.emplace_back("-s");
290         mdArgs_.emplace_back(filename);
291
292         // Set checkpoint file name
293         mdArgs_.emplace_back("-cpi");
294         mdArgs_.emplace_back("state.cpt");
295         /* Note: we normalize the checkpoint file name, but not its full path.
296          * Through version 0.0.8, gmxapi clients change working directory
297          * for each session, so relative path(s) below are appropriate.
298          * A future gmxapi version should avoid changing directories once the
299          * process starts and instead manage files (paths) in an absolute and
300          * immutable way, with abstraction provided through the Context chain-of-responsibility.
301          * TODO: API abstractions for initializing simulations that may be new or partially
302          * complete. Reference gmxapi milestone 13 at https://gitlab.com/gromacs/gromacs/-/issues/2585
303          */
304
305         // Create a mock argv. Note that argv[0] is expected to hold the program name.
306         const int  offset = 1;
307         const auto argc   = static_cast<size_t>(mdArgs_.size() + offset);
308         auto       argv   = std::vector<char*>(argc, nullptr);
309         // argv[0] is ignored, but should be a valid string (e.g. null terminated array of char)
310         argv[0]  = new char[1];
311         *argv[0] = '\0';
312         for (size_t argvIndex = offset; argvIndex < argc; ++argvIndex)
313         {
314             const auto& mdArg = mdArgs_[argvIndex - offset];
315             argv[argvIndex]   = new char[mdArg.length() + 1];
316             strcpy(argv[argvIndex], mdArg.c_str());
317         }
318
319         auto mdModules = std::make_unique<MDModules>();
320
321         const char* desc[] = { "gmxapi placeholder text" };
322
323         // LegacyMdrunOptions needs to be kept alive for the life of ContextImpl,
324         // so we use a data member for now.
325         gmx::LegacyMdrunOptions& options = options_;
326         if (options.updateFromCommandLine(argc, argv.data(), desc) == 0)
327         {
328             return nullptr;
329         }
330
331         ArrayRef<const std::string> multiSimDirectoryNames =
332                 opt2fnsIfOptionSet("-multidir", ssize(options.filenames), options.filenames.data());
333
334
335         // The SimulationContext is necessary with gmxapi so that
336         // resources owned by the client code can have suitable
337         // lifetime. The gmx wrapper binary uses the same infrastructure,
338         // but the lifetime is now trivially that of the invocation of the
339         // wrapper binary.
340         //
341         // For now, this should match the communicator used for hardware
342         // detection. There's no way to assert this is true.
343         auto communicator = mpi_.communicator();
344         // Confirm the precondition for simulationContext().
345         GMX_ASSERT(communicator == MPI_COMM_NULL ? !GMX_LIB_MPI : GMX_LIB_MPI,
346                    "Context communicator does not have an appropriate value for the environment.");
347         SimulationContext simulationContext(communicator, multiSimDirectoryNames);
348
349
350         StartingBehavior startingBehavior        = StartingBehavior::NewSimulation;
351         LogFilePtr       logFileGuard            = nullptr;
352         gmx_multisim_t*  ms                      = simulationContext.multiSimulation_.get();
353         std::tie(startingBehavior, logFileGuard) = handleRestart(
354                 findIsSimulationMasterRank(ms, simulationContext.simulationCommunicator_),
355                 simulationContext.simulationCommunicator_, ms, options.mdrunOptions.appendingBehavior,
356                 ssize(options.filenames), options.filenames.data());
357
358         auto builder = MdrunnerBuilder(std::move(mdModules),
359                                        compat::not_null<SimulationContext*>(&simulationContext));
360         builder.addHardwareDetectionResult(hardwareInformation_.get());
361         builder.addSimulationMethod(options.mdrunOptions, options.pforce, startingBehavior);
362         builder.addDomainDecomposition(options.domdecOptions);
363         // \todo pass by value
364         builder.addNonBonded(options.nbpu_opt_choices[0]);
365         // \todo pass by value
366         builder.addElectrostatics(options.pme_opt_choices[0], options.pme_fft_opt_choices[0]);
367         builder.addBondedTaskAssignment(options.bonded_opt_choices[0]);
368         builder.addUpdateTaskAssignment(options.update_opt_choices[0]);
369         builder.addNeighborList(options.nstlist_cmdline);
370         builder.addReplicaExchange(options.replExParams);
371         // Need to establish run-time values from various inputs to provide a resource handle to Mdrunner
372         builder.addHardwareOptions(options.hw_opt);
373
374         // \todo File names are parameters that should be managed modularly through further factoring.
375         builder.addFilenames(options.filenames);
376         // TODO: Remove `s` and `-cpi` from LegacyMdrunOptions before launch(). #3652
377         auto simulationInput = makeSimulationInput(options);
378         builder.addInput(simulationInput);
379
380         // Note: The gmx_output_env_t life time is not managed after the call to parse_common_args.
381         // \todo Implement lifetime management for gmx_output_env_t.
382         // \todo Output environment should be configured outside of Mdrunner and provided as a resource.
383         builder.addOutputEnvironment(options.oenv);
384         builder.addLogFile(logFileGuard.get());
385
386         // Note, creation is not mature enough to be exposed in the external API yet.
387         launchedSession = createSession(shared_from_this(), std::move(builder),
388                                         std::move(simulationContext), std::move(logFileGuard));
389
390         // Clean up argv once builder is no longer in use
391         // TODO: Remove long-lived references to argv so this is no longer necessary.
392         //       Ref https://gitlab.com/gromacs/gromacs/-/issues/2877
393         for (auto&& string : argv)
394         {
395             if (string != nullptr)
396             {
397                 delete[] string;
398                 string = nullptr;
399             }
400         }
401     }
402     else
403     {
404         throw gmxapi::ProtocolError("Tried to launch a session while a session is still active.");
405     }
406
407     if (launchedSession != nullptr)
408     {
409         // Update weak reference.
410         session_ = launchedSession;
411     }
412     return launchedSession;
413 }
414
415 std::shared_ptr<Session> Context::launch(const Workflow& work)
416 {
417     return impl_->launch(work);
418 }
419
420 Context::Context(std::shared_ptr<ContextImpl> impl) : impl_{ std::move(impl) }
421 {
422     if (!impl_)
423     {
424         throw UsageError("Context requires a non-null implementation member.");
425     }
426 }
427
428 void Context::setMDArgs(const MDArgs& mdArgs)
429 {
430     impl_->mdArgs_ = mdArgs;
431 }
432
433 Context::~Context() = default;
434
435 } // end namespace gmxapi