AmpTools
LikelihoodCalculatorMPI.cc
Go to the documentation of this file.
1 //******************************************************************************
2 // This file is part of AmpTools, a package for performing Amplitude Analysis
3 //
4 // Copyright Trustees of Indiana University 2010, all rights reserved
5 //
6 // This software written by Matthew Shepherd, Ryan Mitchell, and
7 // Hrayr Matevosyan at Indiana University, Bloomington
8 //
9 // Redistribution and use in source and binary forms, with or without
10 // modification, are permitted provided that the following conditions
11 // are met:
12 // 1. Redistributions of source code must retain the above copyright
13 // notice and author attribution, this list of conditions and the
14 // following disclaimer.
15 // 2. Redistributions in binary form must reproduce the above copyright
16 // notice and author attribution, this list of conditions and the
17 // following disclaimer in the documentation and/or other materials
18 // provided with the distribution.
19 // 3. Neither the name of the University nor the names of its contributors
20 // may be used to endorse or promote products derived from this software
21 // without specific prior written permission.
22 //
23 // Creation of derivative forms of this software for commercial
24 // utilization may be subject to restriction; written permission may be
25 // obtained from the Trustees of Indiana University.
26 //
27 // INDIANA UNIVERSITY AND THE AUTHORS MAKE NO REPRESENTATIONS OR WARRANTIES,
28 // EXPRESS OR IMPLIED. By way of example, but not limitation, INDIANA
29 // UNIVERSITY MAKES NO REPRESENTATIONS OR WARRANTIES OF MERCANTABILITY OR
30 // FITNESS FOR ANY PARTICULAR PURPOSE OR THAT THE USE OF THIS SOFTWARE OR
31 // DOCUMENTATION WILL NOT INFRINGE ANY PATENTS, COPYRIGHTS, TRADEMARKS,
32 // OR OTHER RIGHTS. Neither Indiana University nor the authors shall be
33 // held liable for any liability with respect to any claim by the user or
34 // any other party arising from use of the program.
35 //******************************************************************************
36 
37 #include <mpi.h>
38 #include <pthread.h>
39 
42 
45  const NormIntInterface& normInt,
46  DataReader* dataReader,
47  DataReader* bkgReader,
48  ParameterManagerMPI& parManager ) :
49 LikelihoodCalculator( intenManager, normInt, dataReader, bkgReader, parManager ),
50 m_intenManager( intenManager ),
51 m_parManager( parManager ),
52 m_thisId( m_idCounter++ ),
53 m_firstPass( true )
54 {
55  setupMPI();
56 
58 
59  if( !m_isMaster ){
60 
61  // check back in with the master after registration
62  MPI_Send( &m_thisId, 1, MPI_INT, 0, MPITag::kIntSend, MPI_COMM_WORLD );
63  }
64  else{
65 
66  // wait for the workers to check in before proceeding -- if this is
67  // not done, the immediate calls to operator()() could behave
68  // unexpectedly
69 
70  int id;
71  MPI_Status status;
72 
73  for( int i = 1; i < m_numProc; ++i ){
74 
75  MPI_Recv( &id, 1, MPI_INT, i, MPITag::kIntSend, MPI_COMM_WORLD, &status );
76 
77  // ids should match
78  assert( m_thisId == id );
79  }
80  }
81 }
82 
84 
85  // have the master break all of the workers out of their
86  // deliverLikelihood loops
87  if( m_isMaster && m_thisId == kFirstId ){
88 
89  // a two element array to hold commands that go to the workers
90  // the first element is the id of the likelihood calculator
91  // the second element is the flag of the command
92  int cmnd[2];
93  cmnd[0] = m_thisId;
94 
95  // break the likelihood manager out of its loop on the workers
97  for( int i = 1; i < m_numProc; ++i ){
98 
99  MPI_Send( cmnd, 2, MPI_INT, i, MPITag::kIntSend, MPI_COMM_WORLD );
100  }
101  }
102 }
103 
104 double
106 {
107  assert( m_isMaster );
108 
109  MPI_Status status;
110 
111  // a two element array to hold commands that go to the workers
112  // the first element is the id of the likelihood calculator
113  // the second element is the flag of the command
114  int cmnd[2];
115  cmnd[0] = m_thisId;
116 
117  // tell all of the workers to update parameters
118  // note: this is a little inefficient since all instances of the
119  // likelihood calculator share the same parameter manager -- this
120  // cause a little extra overhead in multiple final-state fits
122  for( int i = 1; i < m_numProc; ++i ){
123 
124  MPI_Send( cmnd, 2, MPI_INT, i, MPITag::kIntSend, MPI_COMM_WORLD );
125  }
126 
127  // tell the master to do parameter update
128  m_parManager.updateParameters();
129 
130  // tell all of the workers to send the partial sums
132  for( int i = 1; i < m_numProc; ++i ){
133 
134  MPI_Send( cmnd, 2, MPI_INT, i, MPITag::kIntSend, MPI_COMM_WORLD );
135  }
136 
137  double lnL = 0;
138  double sumBkgWeights = 0;
139  double numBkgEvents = 0;
140  double numDataEvents = 0;
141  double data[4];
142 
143  // collect the sums
144  for( int i = 1; i < m_numProc; ++i ){
145 
146  MPI_Recv( &data, 4, MPI_DOUBLE, i, MPITag::kDoubleSend,
147  MPI_COMM_WORLD, &status );
148 
149  lnL += data[0];
150  sumBkgWeights += data[1];
151  numBkgEvents += data[2];
152  numDataEvents += data[3];
153  }
154 
155  setSumBkgWeights( sumBkgWeights );
156  setNumBkgEvents ( numBkgEvents );
157  setNumDataEvents( numDataEvents );
158 
159  // if we have an amplitude with a free parameter, the call to normIntTerm()
160  // will trigger recomputation of NI's -- we need to put the workers in the
161  // loop to send the recomputed NI's back to the master
162  if( m_intenManager.hasTermWithFreeParam() || m_firstPass ){
163 
165  for( int i = 1; i < m_numProc; ++i ){
166 
167  MPI_Send( cmnd, 2, MPI_INT, i, MPITag::kIntSend, MPI_COMM_WORLD );
168  }
169  }
170 
171  // this call will utilize the NormIntInterface on the master which
172  // ultimately gets handled by the instance of NormIntInterfaceMPI
173  // that is passed into the constructor of this class
174  lnL -= normIntTerm();
175 
176  m_firstPass = false;
177 
178  return -2 * lnL;
179 }
180 
181 void
182 LikelihoodCalculatorMPI::updateParameters()
183 {
184  assert( !m_isMaster );
185 
186  // do the update on the worker nodes
187  m_parManager.updateParameters();
188 }
189 
190 void
191 LikelihoodCalculatorMPI::updateAmpParameter()
192 {
193  assert( !m_isMaster );
194 
195  // do the update on the worker nodes
196  m_parManager.updateAmpParameter();
197 }
198 
199 void
200 LikelihoodCalculatorMPI::computeLikelihood()
201 {
202  assert( !m_isMaster );
203 
204  double data[4];
205 
206  data[0] = dataTerm();
207  data[1] = sumBkgWeights();
208  data[2] = numBkgEvents();
209  data[3] = numDataEvents();
210 
211  MPI_Send( &data, 4, MPI_DOUBLE, 0, MPITag::kDoubleSend, MPI_COMM_WORLD );
212 }
213 
214 void
215 LikelihoodCalculatorMPI::setupMPI()
216 {
217  MPI_Comm_rank( MPI_COMM_WORLD, &m_rank );
218  MPI_Comm_size( MPI_COMM_WORLD, &m_numProc );
219 
220  m_isMaster = ( m_rank == 0 );
221 }
222 
223 int LikelihoodCalculatorMPI::m_idCounter = kFirstId;
void setNumBkgEvents(double num)
void setSumBkgWeights(double sum)
void setNumDataEvents(double num)
static void registerCalculator(int id, LikelihoodCalculatorMPI *calc)
void updateAmpParameter(const string &parName="")
virtual bool hasTermWithFreeParam() const =0
LikelihoodCalculatorMPI(const IntensityManager &intenManager, const NormIntInterface &normInt, DataReader *dataReaderSignal, DataReader *dataReaderBkgnd, ParameterManagerMPI &parManager)
double sumBkgWeights() const