Support ResourceAssignment interface for tMPI library builds.
[gromacs.git] / api / gmxapi / cpp / context.cpp
blob9769d03b84a78fc72c7453a97999c6ad663ce747
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
35 /*! \file
36 * \brief Implementation details of gmxapi::Context
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.
42 * \author M. Eric Irrgang <ericirrgang@gmail.com>
43 * \ingroup gmxapi
46 #include "gmxapi/context.h"
48 #include <cstring>
50 #include <memory>
51 #include <utility>
52 #include <vector>
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"
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"
76 #include "context_impl.h"
77 #include "createsession.h"
78 #include "session_impl.h"
79 #include "workflow.h"
81 namespace gmxapi
84 // Support some tag dispatch. Warning: These are just aliases (not strong types).
85 /*!
86 * \brief Logical helper alias to convert preprocessor constant to type.
88 * \tparam Value Provide the GMX\_LIB\_MPI macro.
90 template<bool Value>
91 using hasLibraryMpi = std::bool_constant<Value>;
92 using gmxThreadMpi = hasLibraryMpi<false>;
93 using gmxLibMpi = hasLibraryMpi<true>;
94 using MpiType = std::conditional_t<GMX_LIB_MPI, gmxLibMpi, gmxThreadMpi>;
96 using MpiContextInitializationError = BasicException<struct MpiContextInitialization>;
99 /*!
100 * \brief Helpers to evaluate the correct precondition for the library build.
102 * TODO: (#3650) Consider distinct MpiContextManager types for clearer definition of preconditions.
104 namespace
107 [[maybe_unused]] MPI_Comm validCommunicator(const MPI_Comm& communicator, const gmxThreadMpi&)
109 if (communicator != MPI_COMM_NULL)
111 throw MpiContextInitializationError(
112 "Provided communicator must be MPI_COMM_NULL for GROMACS built without MPI "
113 "library.");
115 return communicator;
118 [[maybe_unused]] MPI_Comm validCommunicator(const MPI_Comm& communicator, const gmxLibMpi&)
120 if (communicator == MPI_COMM_NULL)
122 throw MpiContextInitializationError("MPI-enabled GROMACS requires a valid communicator.");
124 return communicator;
128 * \brief Return the communicator if it is appropriate for the environment.
130 * \throws MpiContextInitializationError if communicator does not match the
131 * MpiContextManager precondition for the current library configuration.
133 MPI_Comm validCommunicator(const MPI_Comm& communicator)
135 return validCommunicator(communicator, MpiType());
138 //! \brief Provide a reasonable default value.
139 MPI_Comm validCommunicator()
141 return GMX_LIB_MPI ? MPI_COMM_WORLD : MPI_COMM_NULL;
144 } // anonymous namespace
146 MpiContextManager::MpiContextManager(MPI_Comm communicator) :
147 communicator_(std::make_unique<MPI_Comm>(validCommunicator(communicator)))
149 // Safely increments the GROMACS MPI initialization counter after checking
150 // whether the MPI library is already initialized. After this call, MPI_Init
151 // or MPI_Init_thread has been called exactly once.
152 gmx::init(nullptr, nullptr);
153 GMX_RELEASE_ASSERT(!GMX_LIB_MPI || gmx_mpi_initialized(),
154 "MPI should be initialized before reaching this point.");
155 if (this->communicator() != MPI_COMM_NULL)
157 // Synchronise at the point of acquiring a MpiContextManager.
158 gmx_barrier(this->communicator());
162 MpiContextManager::~MpiContextManager()
164 if (communicator_)
166 // This is always safe to call. It is a no-op if
167 // thread-MPI, and if the constructor completed then the
168 // MPI library is initialized with reference counting.
169 gmx::finalize();
173 MpiContextManager::MpiContextManager() : MpiContextManager(validCommunicator()) {}
175 MPI_Comm MpiContextManager::communicator() const
177 if (!communicator_)
179 throw UsageError("Invalid MpiContextManager. Accessed after `move`?");
181 return *communicator_;
184 ContextImpl::~ContextImpl() = default;
186 [[maybe_unused]] static Context createContext(const ResourceAssignment& resources, const gmxLibMpi&)
188 CommHandle handle;
189 resources.applyCommunicator(&handle);
190 if (handle.communicator == MPI_COMM_NULL)
192 throw UsageError("MPI-enabled Simulator contexts require a valid communicator.");
194 auto contextmanager = MpiContextManager(handle.communicator);
195 auto impl = ContextImpl::create(std::move(contextmanager));
196 GMX_ASSERT(impl, "ContextImpl creation method should not be able to return null.");
197 auto context = Context(impl);
198 return context;
201 [[maybe_unused]] static Context createContext(const ResourceAssignment& resources, const gmxThreadMpi&)
203 if (resources.size() > 1)
205 throw UsageError("Only one thread-MPI Simulation per Context is supported.");
207 // Thread-MPI Context does not yet have a need for user-provided resources.
208 // However, see #3650.
209 return createContext();
212 Context createContext(const ResourceAssignment& resources)
214 return createContext(resources, hasLibraryMpi<GMX_LIB_MPI>());
217 Context createContext()
219 MpiContextManager contextmanager;
220 auto impl = ContextImpl::create(std::move(contextmanager));
221 GMX_ASSERT(impl, "ContextImpl creation method should not be able to return null.");
222 auto context = Context(impl);
223 return context;
226 ContextImpl::ContextImpl(MpiContextManager&& mpi) noexcept(std::is_nothrow_constructible_v<gmx::LegacyMdrunOptions>) :
227 mpi_(std::move(mpi))
229 // Confirm our understanding of the MpiContextManager invariant.
230 GMX_ASSERT(mpi_.communicator() == MPI_COMM_NULL ? !GMX_LIB_MPI : GMX_LIB_MPI,
231 "Precondition is an appropriate communicator for the library environment.");
232 // Make sure we didn't change the data members and overlook implementation details.
233 GMX_ASSERT(session_.expired(),
234 "This implementation assumes an expired weak_ptr at initialization.");
237 std::shared_ptr<ContextImpl> ContextImpl::create(MpiContextManager&& mpi)
239 std::shared_ptr<ContextImpl> impl;
240 impl.reset(new ContextImpl(std::move(mpi)));
241 return impl;
244 std::shared_ptr<Session> ContextImpl::launch(const Workflow& work)
246 using namespace gmx;
247 // Much of this implementation is not easily testable: we need tools to inspect simulation
248 // results and to modify simulation inputs.
250 std::shared_ptr<Session> launchedSession = nullptr;
252 // This implementation can only run one workflow at a time.
253 // Check whether we are already aware of an active session.
254 if (session_.expired())
256 // Check workflow spec, build graph for current context, launch and return new session.
257 // \todo This is specific to the session implementation...
258 auto mdNode = work.getNode("MD");
259 std::string filename{};
260 if (mdNode != nullptr)
262 filename = mdNode->params();
265 /* Mock up the argv interface used by option processing infrastructure.
267 * As default behavior, automatically extend trajectories from the checkpoint file.
268 * In the future, our API for objects used to initialize a simulation needs to address the fact that currently a
269 * microstate requires data from both the TPR and checkpoint file to be fully specified. Put another way,
270 * current
271 * GROMACS simulations can take a "configuration" as input that does not constitute a complete microstate in
272 * terms of hidden degrees of freedom (integrator/thermostat/barostat/PRNG state), but we want a clear notion of
273 * a microstate for gmxapi interfaces.
275 * TODO: Remove `-s` and `-cpi` arguments.
276 * Ref: https://gitlab.com/gromacs/gromacs/-/issues/3652
279 // Set input TPR name
280 mdArgs_.emplace_back("-s");
281 mdArgs_.emplace_back(filename);
283 // Set checkpoint file name
284 mdArgs_.emplace_back("-cpi");
285 mdArgs_.emplace_back("state.cpt");
286 /* Note: we normalize the checkpoint file name, but not its full path.
287 * Through version 0.0.8, gmxapi clients change working directory
288 * for each session, so relative path(s) below are appropriate.
289 * A future gmxapi version should avoid changing directories once the
290 * process starts and instead manage files (paths) in an absolute and
291 * immutable way, with abstraction provided through the Context chain-of-responsibility.
292 * TODO: API abstractions for initializing simulations that may be new or partially
293 * complete. Reference gmxapi milestone 13 at https://gitlab.com/gromacs/gromacs/-/issues/2585
296 // Create a mock argv. Note that argv[0] is expected to hold the program name.
297 const int offset = 1;
298 const auto argc = static_cast<size_t>(mdArgs_.size() + offset);
299 auto argv = std::vector<char*>(argc, nullptr);
300 // argv[0] is ignored, but should be a valid string (e.g. null terminated array of char)
301 argv[0] = new char[1];
302 *argv[0] = '\0';
303 for (size_t argvIndex = offset; argvIndex < argc; ++argvIndex)
305 const auto& mdArg = mdArgs_[argvIndex - offset];
306 argv[argvIndex] = new char[mdArg.length() + 1];
307 strcpy(argv[argvIndex], mdArg.c_str());
310 auto mdModules = std::make_unique<MDModules>();
312 const char* desc[] = { "gmxapi placeholder text" };
313 if (options_.updateFromCommandLine(argc, argv.data(), desc) == 0)
315 return nullptr;
318 ArrayRef<const std::string> multiSimDirectoryNames =
319 opt2fnsIfOptionSet("-multidir", ssize(options_.filenames), options_.filenames.data());
321 // The SimulationContext is necessary with gmxapi so that
322 // resources owned by the client code can have suitable
323 // lifetime. The gmx wrapper binary uses the same infrastructure,
324 // but the lifetime is now trivially that of the invocation of the
325 // wrapper binary.
326 auto communicator = mpi_.communicator();
327 // Confirm the precondition for simulationContext().
328 GMX_ASSERT(communicator == MPI_COMM_NULL ? !GMX_LIB_MPI : GMX_LIB_MPI,
329 "Context communicator does not have an appropriate value for the environment.");
330 SimulationContext simulationContext(communicator, multiSimDirectoryNames);
333 StartingBehavior startingBehavior = StartingBehavior::NewSimulation;
334 LogFilePtr logFileGuard = nullptr;
335 gmx_multisim_t* ms = simulationContext.multiSimulation_.get();
336 std::tie(startingBehavior, logFileGuard) =
337 handleRestart(findIsSimulationMasterRank(ms, simulationContext.simulationCommunicator_),
338 communicator, ms, options_.mdrunOptions.appendingBehavior,
339 ssize(options_.filenames), options_.filenames.data());
341 auto builder = MdrunnerBuilder(std::move(mdModules),
342 compat::not_null<SimulationContext*>(&simulationContext));
343 builder.addSimulationMethod(options_.mdrunOptions, options_.pforce, startingBehavior);
344 builder.addDomainDecomposition(options_.domdecOptions);
345 // \todo pass by value
346 builder.addNonBonded(options_.nbpu_opt_choices[0]);
347 // \todo pass by value
348 builder.addElectrostatics(options_.pme_opt_choices[0], options_.pme_fft_opt_choices[0]);
349 builder.addBondedTaskAssignment(options_.bonded_opt_choices[0]);
350 builder.addUpdateTaskAssignment(options_.update_opt_choices[0]);
351 builder.addNeighborList(options_.nstlist_cmdline);
352 builder.addReplicaExchange(options_.replExParams);
353 // Need to establish run-time values from various inputs to provide a resource handle to Mdrunner
354 builder.addHardwareOptions(options_.hw_opt);
356 // \todo File names are parameters that should be managed modularly through further factoring.
357 builder.addFilenames(options_.filenames);
358 // TODO: Remove `s` and `-cpi` from LegacyMdrunOptions before launch(). #3652
359 auto simulationInput = makeSimulationInput(options_);
360 builder.addInput(simulationInput);
362 // Note: The gmx_output_env_t life time is not managed after the call to parse_common_args.
363 // \todo Implement lifetime management for gmx_output_env_t.
364 // \todo Output environment should be configured outside of Mdrunner and provided as a resource.
365 builder.addOutputEnvironment(options_.oenv);
366 builder.addLogFile(logFileGuard.get());
368 // Note, creation is not mature enough to be exposed in the external API yet.
369 launchedSession = createSession(shared_from_this(), std::move(builder),
370 std::move(simulationContext), std::move(logFileGuard));
372 // Clean up argv once builder is no longer in use
373 // TODO: Remove long-lived references to argv so this is no longer necessary.
374 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2877
375 for (auto&& string : argv)
377 if (string != nullptr)
379 delete[] string;
380 string = nullptr;
384 else
386 throw gmxapi::ProtocolError("Tried to launch a session while a session is still active.");
389 if (launchedSession != nullptr)
391 // Update weak reference.
392 session_ = launchedSession;
394 return launchedSession;
397 std::shared_ptr<Session> Context::launch(const Workflow& work)
399 return impl_->launch(work);
402 Context::Context(std::shared_ptr<ContextImpl> impl) : impl_{ std::move(impl) }
404 if (!impl_)
406 throw UsageError("Context requires a non-null implementation member.");
410 void Context::setMDArgs(const MDArgs& mdArgs)
412 impl_->mdArgs_ = mdArgs;
415 Context::~Context() = default;
417 } // end namespace gmxapi