e2efd4f812a88088dc83cdad4a69e923409f41eb
[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/mdlib/stophandler.h"
59 #include "gromacs/mdrunutility/logging.h"
60 #include "gromacs/mdrunutility/multisim.h"
61 #include "gromacs/mdrun/runner.h"
62 #include "gromacs/mdrunutility/handlerestart.h"
63 #include "gromacs/utility/arraysize.h"
64 #include "gromacs/utility/basenetwork.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/gmxmpi.h"
67 #include "gromacs/utility/init.h"
68 #include "gromacs/utility/smalloc.h"
69
70 #include "gmxapi/mpi/multiprocessingresources.h"
71 #include "gmxapi/exceptions.h"
72 #include "gmxapi/session.h"
73 #include "gmxapi/status.h"
74 #include "gmxapi/version.h"
75
76 #include "context_impl.h"
77 #include "createsession.h"
78 #include "session_impl.h"
79 #include "workflow.h"
80
81 namespace gmxapi
82 {
83
84 MpiContextManager::MpiContextManager(MPI_Comm communicator) :
85     communicator_(std::make_unique<MPI_Comm>(communicator))
86 {
87     GMX_ASSERT(*communicator_ == MPI_COMM_NULL ? !GMX_LIB_MPI : GMX_LIB_MPI,
88                "Invalid communicator value for this library configuration.");
89     // Safely increments the GROMACS MPI initialization counter after checking
90     // whether the MPI library is already initialized. After this call, MPI_Init
91     // or MPI_Init_thread has been called exactly once.
92     gmx::init(nullptr, nullptr);
93     GMX_RELEASE_ASSERT(!GMX_LIB_MPI || gmx_mpi_initialized(),
94                        "MPI should be initialized before reaching this point.");
95     if (this->communicator() != MPI_COMM_NULL)
96     {
97         // Synchronise at the point of acquiring a MpiContextManager.
98         gmx_barrier(this->communicator());
99     }
100 };
101
102 MpiContextManager::~MpiContextManager()
103 {
104     if (communicator_)
105     {
106         // This is always safe to call. It is a no-op if
107         // thread-MPI, and if the constructor completed then the
108         // MPI library is initialized with reference counting.
109         gmx::finalize();
110     }
111 }
112
113 MpiContextManager::MpiContextManager() :
114     MpiContextManager(GMX_LIB_MPI ? MPI_COMM_WORLD : MPI_COMM_NULL)
115 {
116 }
117
118 MPI_Comm MpiContextManager::communicator() const
119 {
120     if (!communicator_)
121     {
122         throw UsageError("Invalid MpiContextManager. Accessed after `move`?");
123     }
124     return *communicator_;
125 }
126
127 ContextImpl::~ContextImpl() = default;
128
129 Context createContext(const ResourceAssignment& resources)
130 {
131     CommHandle handle;
132     resources.applyCommunicator(&handle);
133     if (handle.communicator == MPI_COMM_NULL)
134     {
135         throw UsageError("MPI-enabled Simulator contexts require a valid communicator.");
136     }
137     auto contextmanager = MpiContextManager(handle.communicator);
138     auto impl           = ContextImpl::create(std::move(contextmanager));
139     GMX_ASSERT(impl, "ContextImpl creation method should not be able to return null.");
140     auto context = Context(impl);
141     return context;
142 }
143
144 Context createContext()
145 {
146     MpiContextManager contextmanager;
147     auto              impl = ContextImpl::create(std::move(contextmanager));
148     GMX_ASSERT(impl, "ContextImpl creation method should not be able to return null.");
149     auto context = Context(impl);
150     return context;
151 }
152
153 ContextImpl::ContextImpl(MpiContextManager&& mpi) noexcept(std::is_nothrow_constructible_v<gmx::LegacyMdrunOptions>) :
154     mpi_(std::move(mpi))
155 {
156     // Confirm our understanding of the MpiContextManager invariant.
157     GMX_ASSERT(mpi_.communicator() == MPI_COMM_NULL ? !GMX_LIB_MPI : GMX_LIB_MPI,
158                "Precondition is an appropriate communicator for the library environment.");
159     // Make sure we didn't change the data members and overlook implementation details.
160     GMX_ASSERT(session_.expired(),
161                "This implementation assumes an expired weak_ptr at initialization.");
162 }
163
164 std::shared_ptr<ContextImpl> ContextImpl::create(MpiContextManager&& mpi)
165 {
166     std::shared_ptr<ContextImpl> impl;
167     impl.reset(new ContextImpl(std::move(mpi)));
168     return impl;
169 }
170
171 std::shared_ptr<Session> ContextImpl::launch(const Workflow& work)
172 {
173     using namespace gmx;
174     // Much of this implementation is not easily testable: we need tools to inspect simulation
175     // results and to modify simulation inputs.
176
177     std::shared_ptr<Session> launchedSession = nullptr;
178
179     // This implementation can only run one workflow at a time.
180     // Check whether we are already aware of an active session.
181     if (session_.expired())
182     {
183         // Check workflow spec, build graph for current context, launch and return new session.
184         // \todo This is specific to the session implementation...
185         auto        mdNode = work.getNode("MD");
186         std::string filename{};
187         if (mdNode != nullptr)
188         {
189             filename = mdNode->params();
190         }
191
192         /* Mock up the argv interface used by option processing infrastructure.
193          *
194          * As default behavior, automatically extend trajectories from the checkpoint file.
195          * In the future, our API for objects used to initialize a simulation needs to address the fact that currently a
196          * microstate requires data from both the TPR and checkpoint file to be fully specified. Put another way,
197          * current
198          * GROMACS simulations can take a "configuration" as input that does not constitute a complete microstate in
199          * terms of hidden degrees of freedom (integrator/thermostat/barostat/PRNG state), but we want a clear notion of
200          * a microstate for gmxapi interfaces.
201          *
202          * TODO: Remove `-s` and `-cpi` arguments.
203          *       Ref: https://gitlab.com/gromacs/gromacs/-/issues/3652
204          */
205
206         // Set input TPR name
207         mdArgs_.emplace_back("-s");
208         mdArgs_.emplace_back(filename);
209
210         // Set checkpoint file name
211         mdArgs_.emplace_back("-cpi");
212         mdArgs_.emplace_back("state.cpt");
213         /* Note: we normalize the checkpoint file name, but not its full path.
214          * Through version 0.0.8, gmxapi clients change working directory
215          * for each session, so relative path(s) below are appropriate.
216          * A future gmxapi version should avoid changing directories once the
217          * process starts and instead manage files (paths) in an absolute and
218          * immutable way, with abstraction provided through the Context chain-of-responsibility.
219          * TODO: API abstractions for initializing simulations that may be new or partially
220          * complete. Reference gmxapi milestone 13 at https://gitlab.com/gromacs/gromacs/-/issues/2585
221          */
222
223         // Create a mock argv. Note that argv[0] is expected to hold the program name.
224         const int  offset = 1;
225         const auto argc   = static_cast<size_t>(mdArgs_.size() + offset);
226         auto       argv   = std::vector<char*>(argc, nullptr);
227         // argv[0] is ignored, but should be a valid string (e.g. null terminated array of char)
228         argv[0]  = new char[1];
229         *argv[0] = '\0';
230         for (size_t argvIndex = offset; argvIndex < argc; ++argvIndex)
231         {
232             const auto& mdArg = mdArgs_[argvIndex - offset];
233             argv[argvIndex]   = new char[mdArg.length() + 1];
234             strcpy(argv[argvIndex], mdArg.c_str());
235         }
236
237         auto mdModules = std::make_unique<MDModules>();
238
239         const char* desc[] = { "gmxapi placeholder text" };
240         if (options_.updateFromCommandLine(argc, argv.data(), desc) == 0)
241         {
242             return nullptr;
243         }
244
245         ArrayRef<const std::string> multiSimDirectoryNames =
246                 opt2fnsIfOptionSet("-multidir", ssize(options_.filenames), options_.filenames.data());
247
248         // The SimulationContext is necessary with gmxapi so that
249         // resources owned by the client code can have suitable
250         // lifetime. The gmx wrapper binary uses the same infrastructure,
251         // but the lifetime is now trivially that of the invocation of the
252         // wrapper binary.
253         auto communicator = mpi_.communicator();
254         // Confirm the precondition for simulationContext().
255         GMX_ASSERT(communicator == MPI_COMM_NULL ? !GMX_LIB_MPI : GMX_LIB_MPI,
256                    "Context communicator does not have an appropriate value for the environment.");
257         SimulationContext simulationContext(communicator, multiSimDirectoryNames);
258
259
260         StartingBehavior startingBehavior = StartingBehavior::NewSimulation;
261         LogFilePtr       logFileGuard     = nullptr;
262         gmx_multisim_t*  ms               = simulationContext.multiSimulation_.get();
263         std::tie(startingBehavior, logFileGuard) =
264                 handleRestart(findIsSimulationMasterRank(ms, simulationContext.simulationCommunicator_),
265                               communicator, ms, options_.mdrunOptions.appendingBehavior,
266                               ssize(options_.filenames), options_.filenames.data());
267
268         auto builder = MdrunnerBuilder(std::move(mdModules),
269                                        compat::not_null<SimulationContext*>(&simulationContext));
270         builder.addSimulationMethod(options_.mdrunOptions, options_.pforce, startingBehavior);
271         builder.addDomainDecomposition(options_.domdecOptions);
272         // \todo pass by value
273         builder.addNonBonded(options_.nbpu_opt_choices[0]);
274         // \todo pass by value
275         builder.addElectrostatics(options_.pme_opt_choices[0], options_.pme_fft_opt_choices[0]);
276         builder.addBondedTaskAssignment(options_.bonded_opt_choices[0]);
277         builder.addUpdateTaskAssignment(options_.update_opt_choices[0]);
278         builder.addNeighborList(options_.nstlist_cmdline);
279         builder.addReplicaExchange(options_.replExParams);
280         // Need to establish run-time values from various inputs to provide a resource handle to Mdrunner
281         builder.addHardwareOptions(options_.hw_opt);
282
283         // \todo File names are parameters that should be managed modularly through further factoring.
284         builder.addFilenames(options_.filenames);
285         // TODO: Remove `s` and `-cpi` from LegacyMdrunOptions before launch(). #3652
286         auto simulationInput = makeSimulationInput(options_);
287         builder.addInput(simulationInput);
288
289         // Note: The gmx_output_env_t life time is not managed after the call to parse_common_args.
290         // \todo Implement lifetime management for gmx_output_env_t.
291         // \todo Output environment should be configured outside of Mdrunner and provided as a resource.
292         builder.addOutputEnvironment(options_.oenv);
293         builder.addLogFile(logFileGuard.get());
294
295         // Note, creation is not mature enough to be exposed in the external API yet.
296         launchedSession = createSession(shared_from_this(), std::move(builder),
297                                         std::move(simulationContext), std::move(logFileGuard));
298
299         // Clean up argv once builder is no longer in use
300         // TODO: Remove long-lived references to argv so this is no longer necessary.
301         //       Ref https://gitlab.com/gromacs/gromacs/-/issues/2877
302         for (auto&& string : argv)
303         {
304             if (string != nullptr)
305             {
306                 delete[] string;
307                 string = nullptr;
308             }
309         }
310     }
311     else
312     {
313         throw gmxapi::ProtocolError("Tried to launch a session while a session is still active.");
314     }
315
316     if (launchedSession != nullptr)
317     {
318         // Update weak reference.
319         session_ = launchedSession;
320     }
321     return launchedSession;
322 }
323
324 std::shared_ptr<Session> Context::launch(const Workflow& work)
325 {
326     return impl_->launch(work);
327 }
328
329 Context::Context(std::shared_ptr<ContextImpl> impl) : impl_{ std::move(impl) }
330 {
331     if (!impl_)
332     {
333         throw UsageError("Context requires a non-null implementation member.");
334     }
335 }
336
337 void Context::setMDArgs(const MDArgs& mdArgs)
338 {
339     impl_->mdArgs_ = mdArgs;
340 }
341
342 Context::~Context() = default;
343
344 } // end namespace gmxapi