26 #include <queso/InterpolationSurrogateBuilder.h>
29 #include <queso/GslVector.h>
30 #include <queso/GslMatrix.h>
31 #include <queso/InterpolationSurrogateHelper.h>
32 #include <queso/StreamUtilities.h>
39 template<
class V,
class M>
43 m_njobs(this->get_default_data().get_paramDomain().env().numSubEnvironments(), 0)
48 template<
class V,
class M>
52 unsigned int n_values = this->get_default_data().n_values();
53 unsigned int n_workers = this->get_default_data().get_paramDomain().env().numSubEnvironments();
55 unsigned int n_jobs = n_values/n_workers;
56 unsigned int n_leftover = this->get_default_data().n_values() % n_workers;
62 for(
unsigned int n = 0; n < n_workers; n++)
63 this->m_njobs[n] = n_jobs;
68 for(
unsigned int n = 0; n < n_workers; n++)
71 this->m_njobs[n] = n_jobs+1;
73 this->m_njobs[n] = n_jobs;
81 template<
class V,
class M>
84 unsigned int n_begin, n_end;
85 this->set_work_bounds( n_begin, n_end );
88 std::vector<unsigned int> local_n(n_end-n_begin);
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);
96 unsigned int count = 0;
99 V domain_vector(this->get_default_data().get_paramDomain().vectorSpace().zeroVector());
102 std::vector<double> values(this->m_data.size());
104 for(
unsigned int n = n_begin; n < n_end; n++ )
106 this->set_domain_vector( n, domain_vector );
108 this->evaluate_model( domain_vector, values );
112 for(
unsigned int s = 0; s < this->m_data.size(); s++ )
113 local_values[s][count] = values[s];
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) );
125 template<
class V,
class M>
128 unsigned int my_subid = this->get_default_data().get_paramDomain().env().subId();
132 for(
unsigned int n = 0; n < my_subid; n++ )
133 n_begin += m_njobs[n];
135 n_end = n_begin + m_njobs[my_subid];
138 template<
class V,
class M>
140 std::vector<double>& local_values,
146 if( my_subrank == 0 )
148 std::vector<double> all_values(data.
n_values());
150 std::vector<unsigned int> all_indices(data.
n_values());
152 std::vector<int> strides;
153 this->compute_strides( strides );
159 inter0comm.
Gatherv( &local_n[0], local_n.size(), MPI_UNSIGNED,
160 &all_indices[0], &m_njobs[0], &strides[0], MPI_UNSIGNED,
162 "InterpolationSurrogateBuilder::sync_data()",
163 "MpiComm::gatherv() failed!" );
165 inter0comm.
Gatherv( &local_values[0], local_values.size(), MPI_DOUBLE,
166 &all_values[0], &m_njobs[0], &strides[0], MPI_DOUBLE,
168 "InterpolationSurrogateBuilder::sync_data()",
169 "MpiComm::gatherv() failed!" );
178 for(
unsigned int n = 0; n < data.
n_values(); n++ )
179 data.
set_value( all_indices[n], all_values[n] );
187 template<
class V,
class M>
191 std::vector<unsigned int> indices(this->get_default_data().
dim());
195 for(
unsigned int d = 0; d < this->get_default_data().dim(); d++ )
197 domain_vector[d] = this->get_default_data().get_x( d, indices[d] );
201 template<
class V,
class M>
204 unsigned int n_subenvs = this->get_default_data().get_paramDomain().env().numSubEnvironments();
206 strides.resize(n_subenvs);
212 for(
unsigned int n = 1; n < n_subenvs; n++ )
217 stride += this->m_njobs[n-1];
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)
InterpolationSurrogateBuilder()
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)
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.
unsigned int n_values() const
The QUESO MPI Communicator Class.
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.
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's model and populate m_values for the given n_points.
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
void sync_values(unsigned int root)
Sync values across all processors from root processor.
Build interpolation-based surrogate.