BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
MPISimulation.cpp
Go to the documentation of this file.
1 // ************************************************************************** //
2 //
3 // BornAgain: simulate and fit scattering at grazing incidence
4 //
5 //! @file Core/Simulation/MPISimulation.cpp
6 //! @brief Implements class MPISimulation.
7 //!
8 //! @homepage http://www.bornagainproject.org
9 //! @license GNU General Public License v3 or higher (see COPYING)
10 //! @copyright Forschungszentrum Jülich GmbH 2015
11 //! @authors Scientific Computing Group at MLZ Garching
12 //! @authors C. Durniak, M. Ganeva, G. Pospelov, W. Van Herck, J. Wuttke
13 //
14 // ************************************************************************** //
15 
17 
18 #ifdef BORNAGAIN_MPI
19 
21 #include <mpi.h>
22 
23 // -----------------------------------------------------------------------------
24 // MPI support
25 // -----------------------------------------------------------------------------
26 
28 {
29  MPI_Status st;
30 
31  int world_size(0), world_rank(0);
32  MPI_Comm_size(MPI_COMM_WORLD, &world_size);
33  MPI_Comm_rank(MPI_COMM_WORLD, &world_rank);
34 
35  if (world_size == 1) {
36  simulation->runSimulation();
37  return;
38  }
39 
40  SimulationOptions& sim_options = simulation->getOptions();
41  unsigned n_threads = sim_options.getNumberOfThreads();
42  unsigned n_batches = world_size;
43  unsigned current_batch = world_rank;
44  ThreadInfo info;
45  info.n_threads = n_threads;
46  info.n_batches = n_batches;
47  info.current_batch = current_batch;
48  sim_options.setThreadInfo(info);
49  simulation->runSimulation();
50 
51  if (world_rank != 0) {
52  std::vector<double> raw = simulation->rawResults();
53  MPI_Send(&raw[0], raw.size(), MPI_DOUBLE, 0, 0, MPI_COMM_WORLD);
54  }
55  if (world_rank == 0) {
56  auto sum_of_raw = simulation->rawResults();
57  size_t total_size = sum_of_raw.size();
58  for (int i = 1; i < world_size; ++i) {
59  std::vector<double> raw(total_size);
60  MPI_Recv(&raw[0], total_size, MPI_DOUBLE, i, 0, MPI_COMM_WORLD, &st);
61  for (size_t i_raw = 0; i_raw < total_size; ++i_raw)
62  sum_of_raw[i_raw] += raw[i_raw];
63  }
64  simulation->setRawResults(sum_of_raw);
65  }
66 }
67 #else
68 
69 #include <stdexcept>
70 
71 // -----------------------------------------------------------------------------
72 // No MPI support
73 // -----------------------------------------------------------------------------
74 
76 {
77  throw std::runtime_error(
78  "MPISimulation::runSimulation() -> Error! Can't run MPI simulation. "
79  "The package was compiled without MPI support (compile with -DBORNAGAIN_MPI=ON)");
80 }
81 #endif // BORNAGAIN_MPI
Defines class MPISimulation.
Defines class Simulation.
void runSimulation(Simulation *simulation)
Collect the different options for simulation.
unsigned getNumberOfThreads() const
void setThreadInfo(const ThreadInfo &thread_info)
Sets the batch and thread information to be used.
Pure virtual base class of OffSpecularSimulation, GISASSimulation and SpecularSimulation.
Definition: Simulation.h:38
const SimulationOptions & getOptions() const
Definition: Simulation.h:92
void runSimulation()
Run a simulation, possibly averaged over parameter distributions.
Definition: Simulation.cpp:200
virtual std::vector< double > rawResults() const =0
virtual void setRawResults(const std::vector< double > &raw_data)=0
Information to run simulation with dedicated number of threads.
Definition: ThreadInfo.h:21
unsigned current_batch
Definition: ThreadInfo.h:25
unsigned n_batches
Definition: ThreadInfo.h:24
unsigned n_threads
Definition: ThreadInfo.h:23