Move functionality to mdrunutility
[alexxy/gromacs.git] / src / api / cpp / session.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018, 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
36 #include "gmxpre.h"
37
38 #include "config.h"
39
40 #include "gmxapi/session.h"
41
42 #include <memory>
43 #include "gromacs/gmxlib/network.h"
44 #include "gromacs/mdlib/sighandler.h"
45 #include "gromacs/mdrunutility/logging.h"
46 #include "gromacs/mdrunutility/multisim.h"
47 #include "gromacs/restraint/restraintpotential.h"
48 #include "gromacs/utility/gmxassert.h"
49 #include "gromacs/utility/basenetwork.h"
50 #include "gromacs/utility/init.h"
51
52 #include "gmxapi/context.h"
53 #include "gmxapi/exceptions.h"
54 #include "gmxapi/status.h"
55 #include "gmxapi/md/mdmodule.h"
56
57 #include "createsession.h"
58 #include "mdsignals.h"
59 #include "session_impl.h"
60 #include "sessionresources.h"
61
62 namespace gmxapi
63 {
64
65 /*!
66  * \brief Provide RAII management of communications resource state.
67  *
68  * To acquire an MpiContextManager is to have assurance that the GROMACS MPI
69  * environment is ready to use. When the MpiContextManager is released or
70  * goes out of scope, the destructor finalizes the resources.
71  *
72  * \todo Figure out how to manage MPI versus tMPI.
73  * \todo gmx::init should take a subcommunicator rather than use MPI_COMM_WORLD
74  * \todo There is no resource for logging or reporting errors during initialization
75  * \todo Clarify relationship with gmx::SimulationContext.
76  *
77  * \ingroup gmxapi
78  */
79 class MpiContextManager
80 {
81     public:
82         MpiContextManager()
83         {
84             gmx::init(nullptr, nullptr);
85 #ifdef GMX_MPI
86 #if GMX_MPI
87 #if GMX_THREAD_MPI
88             // With thread-MPI, gmx_mpi_initialized() is false until after
89             // spawnThreads in the middle of gmx::Mdrunner::mdrunner(), but
90             // without thread-MPI, MPI is initialized by the time gmx::init()
91             // returns. In other words, this is not an effective context manager
92             // for thread-MPI, but it should be effective for MPI.
93             // \todo Distinguish scope / lifetime for comm resources from implementation details.
94             // \todo Normalize scope / lifetime of comm resources.
95 #else
96             GMX_ASSERT(gmx_mpi_initialized(), "MPI should be initialized before reaching this point.");
97 #endif // GMX_THREAD_MPI
98 #endif // GMX_MPI
99 #endif // defined(GMX_MPI)
100         };
101
102         ~MpiContextManager()
103         {
104             // This is safe to call. It is a no-op if thread-MPI, and if the
105             // constructor completed then MPI is initialized.
106             gmx::finalize();
107         }
108
109         /*!
110          * \brief Exclusive ownership of a scoped context means copying is impossible.
111          *
112          * \{
113          */
114         MpiContextManager(const MpiContextManager &)            = delete;
115         MpiContextManager &operator=(const MpiContextManager &) = delete;
116         //! \}
117
118         /*!
119          * \brief Move semantics are trivial.
120          *
121          * \{
122          */
123         MpiContextManager(MpiContextManager &&) noexcept            = default;
124         MpiContextManager &operator=(MpiContextManager &&) noexcept = default;
125         //! \}
126 };
127
128 SignalManager::SignalManager(gmx::StopHandlerBuilder* stopHandlerBuilder) :
129     state_(std::make_shared<gmx::StopSignal>(gmx::StopSignal::noSignal))
130 {
131
132     /*!
133      * \brief Signal issuer managed by this object.
134      *
135      * Created and bound when the runner is built. Subsequently, client
136      * stop signals are proxied to the simulation through the session
137      * resources. The MD internal signal handler calls this functor
138      * during the MD loop to see if the simulation should be stopped.
139      * Thus, this should execute within a very small fraction of an MD
140      * step and not require any synchronization.
141      */
142     auto currentState     = state_;
143     auto stopSignalIssuer = [currentState](){
144             return *currentState;
145         };
146     stopHandlerBuilder->registerStopCondition(stopSignalIssuer);
147 }
148
149 //! \cond
150 SignalManager::~SignalManager() = default;
151 //! \endcond
152
153 bool SessionImpl::isOpen() const noexcept
154 {
155     // Currently, an active session is equivalent to an active Mdrunner.
156     return bool(runner_);
157 }
158
159 Status SessionImpl::close()
160 {
161     // Assume unsuccessful until proven otherwise.
162     auto successful = Status(false);
163
164     // When the Session is closed, we need to know that the MD output has been finalized.
165     runner_.reset();
166     logFilePtr_.reset();
167
168     // \todo provide meaningful result.
169     // We should be checking that resources were properly shut down, but
170     // there isn't currently a way to do that.
171     successful = true;
172     return successful;
173 }
174
175 Status SessionImpl::run() noexcept
176 {
177     // Status is failure until proven otherwise.
178     auto successful = Status(false);
179     GMX_ASSERT(runner_, "SessionImpl invariant implies valid Mdrunner handle.");
180     auto rc = runner_->mdrunner();
181     if (rc == 0)
182     {
183         successful = true;
184     }
185     return successful;
186 }
187
188 std::unique_ptr<SessionImpl> SessionImpl::create(std::shared_ptr<ContextImpl>  context,
189                                                  gmx::MdrunnerBuilder        &&runnerBuilder,
190                                                  const gmx::SimulationContext &simulationContext,
191                                                  gmx::LogFilePtr               logFilehandle,
192                                                  gmx_multisim_t              * multiSim)
193 {
194     // We should be able to get a communicator (or subcommunicator) through the
195     // Context.
196     return std::make_unique<SessionImpl>(std::move(context),
197                                          std::move(runnerBuilder),
198                                          simulationContext,
199                                          std::move(logFilehandle),
200                                          multiSim);
201 }
202
203 SessionImpl::SessionImpl(std::shared_ptr<ContextImpl>  context,
204                          gmx::MdrunnerBuilder        &&runnerBuilder,
205                          const gmx::SimulationContext &simulationContext,
206                          gmx::LogFilePtr               fplog,
207                          gmx_multisim_t              * multiSim) :
208     context_(std::move(context)),
209     mpiContextManager_(std::make_unique<MpiContextManager>()),
210     simulationContext_(simulationContext),
211     logFilePtr_(std::move(fplog)),
212     multiSim_(multiSim)
213 {
214     GMX_ASSERT(context_, "SessionImpl invariant implies valid ContextImpl handle.");
215     GMX_ASSERT(mpiContextManager_, "SessionImpl invariant implies valid MpiContextManager guard.");
216     GMX_ASSERT(simulationContext_.communicationRecord_, "SessionImpl invariant implies valid commrec.");
217     GMX_UNUSED_VALUE(multiSim_);
218     GMX_UNUSED_VALUE(simulationContext_);
219
220     // \todo Session objects can have logic specialized for the runtime environment.
221
222     auto stopHandlerBuilder = std::make_unique<gmx::StopHandlerBuilder>();
223     signalManager_ = std::make_unique<SignalManager>(stopHandlerBuilder.get());
224     GMX_ASSERT(signalManager_, "SessionImpl invariant includes a valid SignalManager.");
225
226     runnerBuilder.addStopHandlerBuilder(std::move(stopHandlerBuilder));
227     runner_ = std::make_unique<gmx::Mdrunner>(runnerBuilder.build());
228     GMX_ASSERT(runner_, "SessionImpl invariant implies valid Mdrunner handle.");
229
230     // For the libgromacs context, a session should explicitly reset global variables that could
231     // have been set in a previous simulation during the same process.
232     gmx_reset_stop_condition();
233 }
234
235 std::shared_ptr<Session> createSession(std::shared_ptr<ContextImpl>  context,
236                                        gmx::MdrunnerBuilder        &&runnerBuilder,
237                                        const gmx::SimulationContext &simulationContext,
238                                        gmx::LogFilePtr               logFilehandle,
239                                        gmx_multisim_t              * multiSim)
240 {
241     auto newSession = SessionImpl::create(std::move(context),
242                                           std::move(runnerBuilder),
243                                           simulationContext,
244                                           std::move(logFilehandle),
245                                           multiSim);
246     auto launchedSession = std::make_shared<Session>(std::move(newSession));
247     return launchedSession;
248 }
249
250 Status SessionImpl::addRestraint(std::shared_ptr<gmxapi::MDModule> module)
251 {
252     GMX_ASSERT(runner_, "SessionImpl invariant implies valid Mdrunner handle.");
253     Status status {
254         false
255     };
256
257     if (module != nullptr)
258     {
259         const auto &name = module->name();
260         if (restraints_.find(name) == restraints_.end())
261         {
262             auto restraint = module->getRestraint();
263             if (restraint != nullptr)
264             {
265                 restraints_.emplace(std::make_pair(name, restraint));
266                 auto sessionResources = createResources(module);
267                 if (!sessionResources)
268                 {
269                     status = false;
270                 }
271                 else
272                 {
273                     runner_->addPotential(restraint, module->name());
274                     status = true;
275                 }
276             }
277         }
278     }
279     return status;
280 }
281
282
283 SignalManager *SessionImpl::getSignalManager()
284 {
285     SignalManager* ptr = nullptr;
286     if (isOpen())
287     {
288         ptr = signalManager_.get();
289     }
290     return ptr;
291 }
292
293 gmx::Mdrunner *SessionImpl::getRunner()
294 {
295     gmx::Mdrunner * runner = nullptr;
296     if (runner_)
297     {
298         runner = runner_.get();
299     }
300     return runner;
301 }
302
303 gmxapi::SessionResources *SessionImpl::getResources(const std::string &name) const noexcept
304 {
305     gmxapi::SessionResources * resources = nullptr;
306     try
307     {
308         resources = resources_.at(name).get();
309     }
310     catch (const std::out_of_range &e)
311     {
312         // named operation does not have any resources registered.
313     }
314
315     return resources;
316 }
317
318 gmxapi::SessionResources *SessionImpl::createResources(std::shared_ptr<gmxapi::MDModule> module) noexcept
319 {
320     // check if resources already exist for this module
321     // If not, create resources and return handle.
322     // Return nullptr for any failure.
323     gmxapi::SessionResources * resources = nullptr;
324     const auto                &name      = module->name();
325     if (resources_.find(name) == resources_.end())
326     {
327         auto resourcesInstance = std::make_unique<SessionResources>(this, name);
328         resources_.emplace(std::make_pair(name, std::move(resourcesInstance)));
329         resources = resources_.at(name).get();
330         // To do: This should be more dynamic.
331         getSignalManager()->addSignaller(name);
332         if (restraints_.find(name) != restraints_.end())
333         {
334             auto restraintRef = restraints_.at(name);
335             auto restraint    = restraintRef.lock();
336             if (restraint)
337             {
338                 restraint->bindSession(resources);
339             }
340         }
341     }
342     ;
343     return resources;
344 }
345
346 Session::Session(std::unique_ptr<SessionImpl> impl) noexcept
347 {
348     if (impl != nullptr)
349     {
350         impl_ = std::move(impl);
351     }
352     GMX_ASSERT(impl_->isOpen(), "SessionImpl invariant implies valid Mdrunner handle.");
353 }
354
355 Status Session::run() noexcept
356 {
357     GMX_ASSERT(impl_, "Session invariant implies valid implementation object handle.");
358
359     const Status status = impl_->run();
360     return status;
361 }
362
363 Status Session::close()
364 {
365     GMX_ASSERT(impl_, "Session invariant implies valid implementation object handle.");
366
367     auto status = Status(false);
368     if (isOpen())
369     {
370         // \todo catch exceptions when we know what they might be
371         status = impl_->close();
372     }
373
374     return status;
375 }
376
377 Session::~Session()
378 {
379     GMX_ASSERT(impl_, "Session invariant implies valid implementation object handle.");
380     if (isOpen())
381     {
382         try
383         {
384             impl_->close();
385         }
386         catch (const std::exception &)
387         {
388             // \todo find some exception-safe things to do with this via the Context interface.
389         }
390     }
391 }
392
393 bool Session::isOpen() const noexcept
394 {
395     GMX_ASSERT(impl_, "Session invariant implies valid implementation object handle.");
396     const auto result = impl_->isOpen();
397     return result;
398 }
399
400 Status addSessionRestraint(Session                         * session,
401                            std::shared_ptr<gmxapi::MDModule> restraint)
402 {
403     auto status = gmxapi::Status(false);
404
405     if (session != nullptr && restraint != nullptr)
406     {
407         // \todo Improve the external / library API facets
408         // so the public API does not need to offer raw pointers.
409         auto sessionImpl = session->getRaw();
410
411         GMX_RELEASE_ASSERT(sessionImpl, "Session invariant implies valid implementation object handle.");
412         // GMX_ASSERT alone is not strong enough to convince linters not to warn of possible nullptr.
413         if (sessionImpl)
414         {
415             status = sessionImpl->addRestraint(std::move(restraint));
416         }
417     }
418     return status;
419 }
420
421 //! \cond internal
422 SessionImpl *Session::getRaw() const noexcept
423 {
424     return impl_.get();
425 }
426 //! \endcond
427
428 std::shared_ptr<Session> launchSession(Context* context, const Workflow &work) noexcept
429 {
430     auto session = context->launch(work);
431     return session;
432 }
433
434 SessionImpl::~SessionImpl() = default;
435
436 SessionResources::SessionResources(gmxapi::SessionImpl *session,
437                                    std::string          name) :
438     sessionImpl_(session),
439     name_(std::move(name))
440 {
441 }
442
443 SessionResources::~SessionResources() = default;
444
445 const std::string SessionResources::name() const
446 {
447     return name_;
448 }
449
450 Signal SessionResources::getMdrunnerSignal(md::signals signal)
451 {
452     //// while there is only one choice...
453     if (signal != md::signals::STOP)
454     {
455         throw gmxapi::NotImplementedError("This signaller only handles stop signals.");
456     }
457
458     // Get a signalling proxy for the caller.
459     auto signalManager = sessionImpl_->getSignalManager();
460     if (signalManager == nullptr)
461     {
462         throw gmxapi::ProtocolError("Client requested access to a signaller that is not available.");
463     }
464     auto functor = signalManager->getSignal(name_, signal);
465
466     return functor;
467 }
468
469 } // end namespace gmxapi