queso-0.53.0
InterpolationSurrogateBuilder.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // QUESO - a library to support the Quantification of Uncertainty
5 // for Estimation, Simulation and Optimization
6 //
7 // Copyright (C) 2008-2015 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 // This class
26 #include <queso/InterpolationSurrogateBuilder.h>
27 
28 // QUESO
29 #include <queso/GslVector.h>
30 #include <queso/GslMatrix.h>
31 #include <queso/InterpolationSurrogateHelper.h>
32 #include <queso/StreamUtilities.h>
33 
34 // C++
35 #include <numeric>
36 
37 namespace QUESO
38 {
39  template<class V, class M>
41  : SurrogateBuilderBase<V>(),
42  m_data(data),
43  m_njobs(this->get_default_data().get_paramDomain().env().numSubEnvironments(), 0)
44  {
45  this->partition_work();
46  }
47 
48  template<class V, class M>
50  {
51  // Convenience
52  unsigned int n_values = this->get_default_data().n_values();
53  unsigned int n_workers = this->get_default_data().get_paramDomain().env().numSubEnvironments();
54 
55  unsigned int n_jobs = n_values/n_workers;
56  unsigned int n_leftover = this->get_default_data().n_values() % n_workers;
57 
58  /* If the number of values is evenly divisible over all workers,
59  then everyone gets the same amount work */
60  if( n_leftover == 0 )
61  {
62  for(unsigned int n = 0; n < n_workers; n++)
63  this->m_njobs[n] = n_jobs;
64  }
65  /* Otherwise, some workers get more work than others*/
66  else
67  {
68  for(unsigned int n = 0; n < n_workers; n++)
69  {
70  if( n < n_leftover )
71  this->m_njobs[n] = n_jobs+1;
72  else
73  this->m_njobs[n] = n_jobs;
74  }
75  }
76 
77  // Sanity check
78  queso_assert_equal_to( n_values, std::accumulate( m_njobs.begin(), m_njobs.end(), 0 ) );
79  }
80 
81  template<class V, class M>
83  {
84  unsigned int n_begin, n_end;
85  this->set_work_bounds( n_begin, n_end );
86 
87  // Cache each processors work, then we only need to do 1 Allgather
88  std::vector<unsigned int> local_n(n_end-n_begin);
89 
90  // We need to cache (n_end-n_begin) values for each dataset,
91  std::vector<std::vector<double> > local_values(this->m_data.size());
92  for( std::vector<std::vector<double> >::iterator it = local_values.begin();
93  it != local_values.end(); ++it )
94  it->resize(n_end-n_begin);
95 
96  unsigned int count = 0;
97 
98  // vector to store current domain value
99  V domain_vector(this->get_default_data().get_paramDomain().vectorSpace().zeroVector());
100 
101  // vector to store values evaluated at the current domain_vector
102  std::vector<double> values(this->m_data.size());
103 
104  for( unsigned int n = n_begin; n < n_end; n++ )
105  {
106  this->set_domain_vector( n, domain_vector );
107 
108  this->evaluate_model( domain_vector, values );
109 
110  local_n[count] = n;
111 
112  for( unsigned int s = 0; s < this->m_data.size(); s++ )
113  local_values[s][count] = values[s];
114 
115  count += 1;
116  }
117 
118  /* Sync all the locally computed values between the subenvironments
119  so all processes have all the computed values. We need to sync
120  values for every data set. */
121  for( unsigned int s = 0; s < this->m_data.size(); s++ )
122  this->sync_data( local_n, local_values[s], this->m_data.get_dataset(s) );
123  }
124 
125  template<class V, class M>
126  void InterpolationSurrogateBuilder<V,M>::set_work_bounds( unsigned int& n_begin, unsigned int& n_end ) const
127  {
128  unsigned int my_subid = this->get_default_data().get_paramDomain().env().subId();
129 
130  /* Starting index will be the sum of the all the previous num jobs */
131  n_begin = 0;
132  for( unsigned int n = 0; n < my_subid; n++ )
133  n_begin += m_njobs[n];
134 
135  n_end = n_begin + m_njobs[my_subid];
136  }
137 
138  template<class V, class M>
139  void InterpolationSurrogateBuilder<V,M>::sync_data( std::vector<unsigned int>& local_n,
140  std::vector<double>& local_values,
142  {
143  // Only members of the inter0comm will do the communication of the local values
144  unsigned int my_subrank = data.get_paramDomain().env().subRank();
145 
146  if( my_subrank == 0 )
147  {
148  std::vector<double> all_values(data.n_values());
149 
150  std::vector<unsigned int> all_indices(data.n_values());
151 
152  std::vector<int> strides;
153  this->compute_strides( strides );
154 
155  const MpiComm& inter0comm = data.get_paramDomain().env().inter0Comm();
156 
159  inter0comm.Gatherv( &local_n[0], local_n.size(), MPI_UNSIGNED,
160  &all_indices[0], &m_njobs[0], &strides[0], MPI_UNSIGNED,
161  0 /*root*/,
162  "InterpolationSurrogateBuilder::sync_data()",
163  "MpiComm::gatherv() failed!" );
164 
165  inter0comm.Gatherv( &local_values[0], local_values.size(), MPI_DOUBLE,
166  &all_values[0], &m_njobs[0], &strides[0], MPI_DOUBLE,
167  0 /*root*/,
168  "InterpolationSurrogateBuilder::sync_data()",
169  "MpiComm::gatherv() failed!" );
170 
171  // Now set the values.
172  /* PB: Although we are guaranteed per-rank ordering of the data we gathered,
173  I'm not sure we can assume the same continuity of the inter0 ranks, i.e.
174  I'm not sure how QUESO ordered the inter0 ranks. So, we go ahead and
175  manually set the values. */
176  if( data.get_paramDomain().env().subRank() == 0 )
177  {
178  for( unsigned int n = 0; n < data.n_values(); n++ )
179  data.set_value( all_indices[n], all_values[n] );
180  }
181  }
182 
183  // Now broadcast the values data to all other processes
184  data.sync_values( 0 /*root*/);
185  }
186 
187  template<class V, class M>
188  void InterpolationSurrogateBuilder<V,M>::set_domain_vector( unsigned int n, V& domain_vector ) const
189  {
190  // Convert global index n to local coordinates in each dimension
191  std::vector<unsigned int> indices(this->get_default_data().dim());
192  InterpolationSurrogateHelper::globalToCoord( n, this->get_default_data().get_n_points(), indices );
193 
194  // Use indices to get x coordinates and populate domain_vector
195  for( unsigned int d = 0; d < this->get_default_data().dim(); d++ )
196  {
197  domain_vector[d] = this->get_default_data().get_x( d, indices[d] );
198  }
199  }
200 
201  template<class V, class M>
202  void InterpolationSurrogateBuilder<V,M>::compute_strides( std::vector<int>& strides ) const
203  {
204  unsigned int n_subenvs = this->get_default_data().get_paramDomain().env().numSubEnvironments();
205 
206  strides.resize(n_subenvs);
207 
208  // Don't stride the first entry
209  strides[0] = 0;
210 
211  int stride = 0;
212  for( unsigned int n = 1; n < n_subenvs; n++ )
213  {
214  // The stride is measured agaisnt the beginning of the buffer
215  // We want things packed tightly together so just stride
216  // by the number of entries from the previous group.
217  stride += this->m_njobs[n-1];
218  strides[n] = stride;
219  }
220  }
221 
222 } // end namespace QUESO
223 
224 // Instantiate
const BoxSubset< V, M > & get_paramDomain() const
static void globalToCoord(unsigned int global, const std::vector< unsigned int > &n_points, std::vector< unsigned int > &coord_indices)
Inverse of coordToGlobal map.
void set_value(unsigned int n, double value)
void partition_work()
Partition the workload of model evaluations across the subenvironments.
void set_domain_vector(unsigned int n, V &domain_vector) const
Provide the spatial coordinates for the global index n.
#define queso_assert_equal_to(expr1, expr2)
Definition: asserts.h:149
void sync_data(std::vector< unsigned int > &local_n, std::vector< double > &local_values, InterpolationSurrogateData< V, M > &data)
Take the local values computed from each process and communicate.
The QUESO MPI Communicator Class.
Definition: MpiComm.h:75
Container class for multiple, consistent InterpolationSurrogateData objects.
void Gatherv(void *sendbuf, int sendcnt, RawType_MPI_Datatype sendtype, void *recvbuf, int *recvcnts, int *displs, RawType_MPI_Datatype recvtype, int root, const char *whereMsg, const char *whatMsg) const
Gathers into specified locations from all processes in a group.
Definition: MpiComm.C:158
void set_work_bounds(unsigned int &n_begin, unsigned int &n_end) const
Set the starting and ending global indices for the current subenvironment.
Base class for builders of surrogates.
void build_values()
Execute the user&#39;s model and populate m_values for the given n_points.
int dim
Definition: ann2fig.cpp:81
void compute_strides(std::vector< int > &strides) const
Helper function to compute strides needed for MPI_Gatherv.
that you receive source code or can get it if you want it
Definition: License.txt:41
void sync_values(unsigned int root)
Sync values across all processors from root processor.
Build interpolation-based surrogate.

Generated on Thu Jun 11 2015 13:52:33 for queso-0.53.0 by  doxygen 1.8.5