Doxygen
1.9.1
|
Typedefs | |
using | TasGrid::ModelSignatureMPI = std::function< void(std::vector< double > const &x, std::vector< double > &y)> |
Signature of a model function to be used in the MPI construction procedures. More... | |
Functions | |
template<bool use_initial_guess = no_initial_guess> | |
void | TasGrid::mpiConstructSurrogate (ModelSignatureMPI model, int num_dimensions, int num_outputs, size_t max_num_points, size_t max_samples_per_job, size_t max_num_ranks, int tagx, int tagy, int root, MPI_Comm comm, TasmanianSparseGrid &grid, double tolerance, TypeRefinement criteria, int output=-1, std::vector< int > const &level_limits=std::vector< int >(), std::vector< double > const &scale_correction=std::vector< double >(), std::string const &checkpoint_filename=std::string()) |
Construct a sparse grid surrogate to the model defined by the lambda, MPI version. More... | |
template<bool use_initial_guess = no_initial_guess> | |
void | TasGrid::mpiConstructSurrogate (ModelSignatureMPI model, int num_dimensions, int num_outputs, size_t max_num_points, size_t max_samples_per_job, size_t max_num_ranks, int tagx, int tagy, int root, MPI_Comm comm, TasmanianSparseGrid &grid, TypeDepth type, std::vector< int > const &anisotropic_weights=std::vector< int >(), std::vector< int > const &level_limits=std::vector< int >(), std::string const &checkpoint_filename=std::string()) |
Construct a sparse grid surrogate to the model defined by the lambda, MPI version. More... | |
template<bool use_initial_guess = no_initial_guess> | |
void | TasGrid::mpiConstructSurrogate (ModelSignatureMPI model, int num_dimensions, int num_outputs, size_t max_num_points, size_t max_samples_per_job, size_t max_num_ranks, int tagx, int tagy, int root, MPI_Comm comm, TasmanianSparseGrid &grid, TypeDepth type, int output, std::vector< int > const &level_limits=std::vector< int >(), std::string const &checkpoint_filename=std::string()) |
Construct a sparse grid surrogate to the model defined by the lambda, MPI version. More... | |
template<bool use_initial_guess = no_initial_guess> | |
void | TasGrid::mpiConstructWorker (ModelSignatureMPI model, int num_dimensions, int num_outputs, size_t max_samples_per_job, size_t max_num_ranks, int tagx, int tagy, int root, MPI_Comm comm) |
Executes the worker (i.e., non-root) portion of the MPI sampling. More... | |
Distributed sampling procedure, construction similar to TasGrid::constructSurrogate() but the parallelism is associated with different ranks in an MPI communicator. These templates require Tasmanian_ENABLE_MPI=ON.
using TasGrid::ModelSignatureMPI = typedef std::function<void(std::vector<double> const &x, std::vector<double> &y)> |
Signature of a model function to be used in the MPI construction procedures.
The logic is identical to the TasGrid::ModelSignature used in the TasGrid::constructSurrogate() templates, with the exception of the removal of the thread-id. In an MPI context, the thread-id is replaced by the MPI rank which can be queried using MPI_Comm_rank().
void TasGrid::mpiConstructSurrogate | ( | ModelSignatureMPI | model, |
int | num_dimensions, | ||
int | num_outputs, | ||
size_t | max_num_points, | ||
size_t | max_samples_per_job, | ||
size_t | max_num_ranks, | ||
int | tagx, | ||
int | tagy, | ||
int | root, | ||
MPI_Comm | comm, | ||
TasmanianSparseGrid & | grid, | ||
double | tolerance, | ||
TypeRefinement | criteria, | ||
int | output = -1 , |
||
std::vector< int > const & | level_limits = std::vector<int>() , |
||
std::vector< double > const & | scale_correction = std::vector<double>() , |
||
std::string const & | checkpoint_filename = std::string() |
||
) |
Construct a sparse grid surrogate to the model defined by the lambda, MPI version.
The logic is the same as in TasGrid::constructSurrogate(), construct adaptive sparse grid surrogate to the mode defined by the lambda using parallel computations. The model will be executed with different inputs x on different ranks within the MPI communicator. All ranks of across the communicator must make an identical call with identical parameters before the grid. The grid and the refinement and checkpoint parameters will be used only by the root rank and will not be addressed by the rest.
The list of inputs is very similar to TasGrid::constructSurrogate(), here we outline the main differences.
use_initial_guess | defines whether to compute an initial guess for the model outputs and send that over the communicator. The MPI procedure always runs in parallel mode hence there is no additional template parameter. |
model | the inputs x and outputs y of the model are the same as in TasGrid::constructSurrogate(); however, there is no thread_id. Only one call of the model will be executed on each active rank and thus the MPI native MPI_Comm_rank() can be used in place of the thread_id. |
num_dimensions | must match grid.getNumDimensions(). This is a required parameter since only the root rank is assumed to have access to the grid, the other ranks need to know the correct dimension without the grid. |
num_outputs | must match grid.getNumOutputs(), see num_dimensions. |
max_num_points | defines the computational budget, same as TasGrid::constructSurrogate(). |
max_samples_per_job | defines the number of samples that will be given to the model in each call, see TasGrid::constructSurrogate(). |
max_num_ranks | defines the number of ranks to use, if it exceeds the size of the communicator then all ranks will be used in the process. If less than the size, the root process will not make calls to the model but it will be working on the internal grid methods. This is useful when dealing with many points where the cost of loading points is not negligible. If the max_num_ranks is less than the communicator size minus one, then some of the ranks will be completely idle (usually the ranks with higher number). |
tagx | is the MPI_Send() tag to use when communicating the inputs x to the worker rank. |
tagy | is the MPI_Send() tag to use when communicating the outputs y to the root rank. |
root | is the MPI rank within the communicator that has the initial copy of the grid as well as the refinement parameters that are identical to TasGrid::constructSurrogate(). The final result from the construction will be held by the root rank only, the grids for the rest of other ranks will not be accessed. |
comm | is the MPI communicator where sampling will take place, e.g., MPI_COMM_WORLD. |
grid | on rank root is the grid to be constructed, ignored when the rank is different. |
tolerance | refinement parameter, see TasGrid::constructSurrogate(). |
criteria | refinement parameter, see TasGrid::constructSurrogate(). |
output | refinement parameter, see TasGrid::constructSurrogate(). |
level_limits | refinement parameter, see TasGrid::constructSurrogate(). |
scale_correction | refinement parameter, see TasGrid::constructSurrogate(). |
checkpoint_filename | is either empty (skip checkpoints) or the filename to use to save the intermediate results. Note that only the root rank will issue read/write commands to the files. |
void TasGrid::mpiConstructSurrogate | ( | ModelSignatureMPI | model, |
int | num_dimensions, | ||
int | num_outputs, | ||
size_t | max_num_points, | ||
size_t | max_samples_per_job, | ||
size_t | max_num_ranks, | ||
int | tagx, | ||
int | tagy, | ||
int | root, | ||
MPI_Comm | comm, | ||
TasmanianSparseGrid & | grid, | ||
TypeDepth | type, | ||
std::vector< int > const & | anisotropic_weights = std::vector<int>() , |
||
std::vector< int > const & | level_limits = std::vector<int>() , |
||
std::string const & | checkpoint_filename = std::string() |
||
) |
Construct a sparse grid surrogate to the model defined by the lambda, MPI version.
Uses the user provided anisotropic_weights to order the samples by importance, see
TasmanianSparseGrid::getCandidateConstructionPoints(). and the overloads to TasGrid::constructSurrogate().
void TasGrid::mpiConstructSurrogate | ( | ModelSignatureMPI | model, |
int | num_dimensions, | ||
int | num_outputs, | ||
size_t | max_num_points, | ||
size_t | max_samples_per_job, | ||
size_t | max_num_ranks, | ||
int | tagx, | ||
int | tagy, | ||
int | root, | ||
MPI_Comm | comm, | ||
TasmanianSparseGrid & | grid, | ||
TypeDepth | type, | ||
int | output, | ||
std::vector< int > const & | level_limits = std::vector<int>() , |
||
std::string const & | checkpoint_filename = std::string() |
||
) |
Construct a sparse grid surrogate to the model defined by the lambda, MPI version.
Uses anisotropic weights to order the samples by importance, starts with a fully isotropic grid until enough points are loaded to allow to estimate the weights, see TasmanianSparseGrid::getCandidateConstructionPoints(). and the overloads to TasGrid::constructSurrogate().
void TasGrid::mpiConstructWorker | ( | ModelSignatureMPI | model, |
int | num_dimensions, | ||
int | num_outputs, | ||
size_t | max_samples_per_job, | ||
size_t | max_num_ranks, | ||
int | tagx, | ||
int | tagy, | ||
int | root, | ||
MPI_Comm | comm | ||
) |
Executes the worker (i.e., non-root) portion of the MPI sampling.
Using any of the TasGrid::mpiConstructSurrogate() overloads, the grid, the budget and the refinement parameters are accessed only by the root rank, the rest is not used by the workers. This template can be instantiated for the non-root rank with only the relevant inputs.