1 /* _______________________________________________________________________
2
3 DAKOTA: Design Analysis Kit for Optimization and Terascale Applications
4 Copyright 2014-2020 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
5 This software is distributed under the GNU Lesser General Public License.
6 For more information, see the README file in the top Dakota directory.
7 _______________________________________________________________________ */
8
9 //- Class: DataFitSurrModel
10 //- Description: A model which manages a surrogate relationship between a
11 //- data fit approximation (local, multipoint, or global) and
12 //- a truth model.
13 //- Owner: Mike Eldred
14 //- Checked by:
15 //- Version: $Id: DataFitSurrModel.hpp 7029 2010-10-22 00:17:02Z mseldre $
16
17 #ifndef DATA_FIT_SURR_MODEL_H
18 #define DATA_FIT_SURR_MODEL_H
19
20 #include "dakota_data_types.hpp"
21 #include "SurrogateData.hpp"
22 #include "SurrogateModel.hpp"
23 #include "DakotaInterface.hpp"
24 #include "DakotaIterator.hpp"
25 #include "ProblemDescDB.hpp"
26 #include "ParallelLibrary.hpp"
27
28
29 namespace Dakota {
30
31 /// Derived model class within the surrogate model branch for managing
32 /// data fit surrogates (global and local)
33
34 /** The DataFitSurrModel class manages global or local approximations
35 (surrogates that involve data fits) that are used in place of an
36 expensive model. The class contains an approxInterface (required
37 for both global and local) which manages the approximate function
38 evaluations, an actualModel (optional for global, required for
39 local) which provides truth evaluations for building the
40 surrogate, and a daceIterator (optional for global, not used for
41 local) which selects parameter sets on which to evaluate
42 actualModel in order to generate the necessary data for building
43 global approximations. */
44
45 class DataFitSurrModel: public SurrogateModel
46 {
47 public:
48
49 //
50 //- Heading: Constructors and destructor
51 //
52
53 /// constructor
54 DataFitSurrModel(ProblemDescDB& problem_db);
55 /// alternate constructor for instantiations on the fly
56 DataFitSurrModel(Iterator& dace_iterator, Model& actual_model,
57 const ActiveSet& set, const String& approx_type,
58 const UShortArray& approx_order, short corr_type,
59 short corr_order, short data_order, short output_level,
60 const String& point_reuse,
61 const String& import_build_points_file = String(),
62 unsigned short import_build_format = TABULAR_ANNOTATED,
63 bool import_build_active_only = false,
64 const String& export_approx_points_file = String(),
65 unsigned short export_approx_format = TABULAR_ANNOTATED);
66 /// destructor
67 ~DataFitSurrModel();
68
69 //
70 //- Heading: Member functions
71 //
72
73 /// set pointsTotal and pointsManagement mode
74 void total_points(int points);
75 /// return points required for build according to pointsManagement mode
76 int required_points();
77
78 void declare_sources();
79
80 protected:
81
82 //
83 //- Heading: Virtual function redefinitions
84 //
85
86 size_t qoi() const;
87
88 DiscrepancyCorrection& discrepancy_correction();
89 short correction_type();
90 void correction_type(short corr_type);
91
92 bool initialize_mapping(ParLevLIter pl_iter);
93 bool finalize_mapping();
94
95 void nested_variable_mappings(const SizetArray& c_index1,
96 const SizetArray& di_index1,
97 const SizetArray& ds_index1,
98 const SizetArray& dr_index1,
99 const ShortArray& c_target2,
100 const ShortArray& di_target2,
101 const ShortArray& ds_target2,
102 const ShortArray& dr_target2);
103 const SizetArray& nested_acv1_indices() const;
104 const ShortArray& nested_acv2_targets() const;
105 short query_distribution_parameter_derivatives() const;
106
107 void check_submodel_compatibility(const Model& sub_model);
108
109 // Perform the response computation portions specific to this derived
110 // class. In this case, it simply employs approxInterface.map()/synch()/
111 // synch_nowait() where approxInterface is a local, multipoint, or global
112 // approximation.
113
114 void derived_evaluate(const ActiveSet& set);
115 void derived_evaluate_nowait(const ActiveSet& set);
116 const IntResponseMap& derived_synchronize();
117 const IntResponseMap& derived_synchronize_nowait();
118
119 /// map incoming ASV into actual request for surrogate construction, managing
120 /// any mismatch in sizes due to response aggregation modes in actualModel
121 void asv_inflate_build(const ShortArray& orig_asv, ShortArray& actual_asv);
122 /// split incoming ASV into actual and approximate evaluation requests,
123 /// managing any mismatch in sizes due to response aggregation modes in
124 /// actualModel
125 void asv_split_eval(const ShortArray& orig_asv, ShortArray& actual_asv,
126 ShortArray& approx_asv);
127
128 /// return daceIterator
129 Iterator& subordinate_iterator();
130
131 /// set active model key within approxInterface
132 void active_model_key(const UShortArray& mi_key);
133 /// remove all model keys within approxInterface
134 void clear_model_keys();
135
136 /// return this model instance
137 Model& surrogate_model();
138 /// return this model instance
139 const Model& surrogate_model() const;
140 /// return actualModel
141 Model& truth_model();
142 /// return actualModel
143 const Model& truth_model() const;
144
145 /// return actualModel (and optionally its sub-models)
146 void derived_subordinate_models(ModelList& ml, bool recurse_flag);
147 /// pass request to actualModel if recursing
148 void resize_from_subordinate_model(size_t depth =
149 std::numeric_limits<size_t>::max());
150 /// pass request to actualModel if recursing and then update from it
151 void update_from_subordinate_model(size_t depth =
152 std::numeric_limits<size_t>::max());
153 /// return approxInterface
154 Interface& derived_interface();
155
156 /// set the relative weightings for multiple objective functions or least
157 /// squares terms and optionally recurses into actualModel
158 void primary_response_fn_weights(const RealVector& wts,
159 bool recurse_flag = true);
160
161 /// set responseMode and pass any bypass request on to actualModel for
162 /// any lower-level surrogates.
163 void surrogate_response_mode(short mode);
164
165 // link together more than one SurrogateData instance within
166 // approxInterface.functionSurfaces[i].approxData[j]
167 //void link_multilevel_approximation_data();
168
169 /// (re)set the surrogate index set in SurrogateModel::surrogateFnIndices
170 /// and ApproximationInterface::approxFnIndices
171 void surrogate_function_indices(const IntSet& surr_fn_indices);
172
173 /// Builds the local/multipoint/global approximation using
174 /// daceIterator/actualModel to generate new data points
175 void build_approximation();
176 /// Builds the local/multipoint/global approximation using
177 /// daceIterator/actualModel to generate new data points that
178 /// augment the passed vars/response anchor point
179 bool build_approximation(const Variables& vars,
180 const IntResponsePair& response_pr);
181
182 /// Rebuilds the local/multipoint/global approximation using
183 /// daceIterator/actualModel to generate an increment of appended data
184 void rebuild_approximation();
185 /// Rebuilds the local/multipoint/global approximation using
186 /// the passed response data for a single sample
187 void rebuild_approximation(const IntResponsePair& response_pr);
188 /// Rebuilds the local/multipoint/global approximation using
189 /// the passed response data for a set of samples
190 void rebuild_approximation(const IntResponseMap& all_resp);
191
192 /// replaces the approximation data with daceIterator results and
193 /// rebuilds the approximation if requested
194 void update_approximation(bool rebuild_flag);
195 /// replaces the anchor point, and rebuilds the approximation if requested
196 void update_approximation(const Variables& vars,
197 const IntResponsePair& response_pr,
198 bool rebuild_flag);
199 /// replaces the current points array and rebuilds the approximation
200 /// if requested
201 void update_approximation(const VariablesArray& vars_array,
202 const IntResponseMap& resp_map, bool rebuild_flag);
203 /// replaces the current points array and rebuilds the approximation
204 /// if requested
205 void update_approximation(const RealMatrix& samples,
206 const IntResponseMap& resp_map, bool rebuild_flag);
207
208 /// appends daceIterator results to a global approximation and rebuilds
209 /// it if requested
210 void append_approximation(bool rebuild_flag);
211 /// appends a point to a global approximation and rebuilds it if requested
212 void append_approximation(const Variables& vars,
213 const IntResponsePair& response_pr,
214 bool rebuild_flag);
215 /// appends an array of points to a global approximation and rebuilds it
216 /// if requested
217 void append_approximation(const VariablesArray& vars_array,
218 const IntResponseMap& resp_map, bool rebuild_flag);
219 /// appends a matrix of points to a global approximation and rebuilds it
220 /// if requested
221 void append_approximation(const RealMatrix& samples,
222 const IntResponseMap& resp_map, bool rebuild_flag);
223
224 /// remove approximation data added on previous append_approximation() call
225 /// or a specified number of points
226 void pop_approximation(bool save_surr_data, bool rebuild_flag = false);
227
228 /// retrieve a previous approximation data state
229 void push_approximation();
230 /// query for whether a trial increment can be retrieved
231 bool push_available();
232
233 /// finalize data fit by applying all previous trial increments
234 void finalize_approximation();
235
236 /// combine all level approximations into a separate composite approximation
237 void combine_approximation();
238 /// promote the combined approximation into the active one
239 void combined_to_active(bool clear_combined = true);
240
241 /// clear inactive data stored in the approxInterface
242 void clear_inactive();
243
244 /// query approxInterface for available advancements in order, rank, etc.
245 bool advancement_available();
246 /// query approxInterface for updates in formulation (requiring a rebuild)
247 bool formulation_updated() const;
248 /// update the formulation status in approxInterface
249 void formulation_updated(bool update);
250
251 /// execute the DACE iterator to generate build data
252 void run_dace();
253
254 /// retrieve the SharedApproxData from approxInterface
255 SharedApproxData& shared_approximation();
256 /// retrieve the set of Approximations from approxInterface
257 std::vector<Approximation>& approximations();
258 /// return the approximation coefficients from each Approximation
259 /// (request forwarded to approxInterface)
260 const RealVectorArray& approximation_coefficients(bool normalized = false);
261 /// set the approximation coefficients within each Approximation
262 /// (request forwarded to approxInterface)
263 void approximation_coefficients(const RealVectorArray& approx_coeffs,
264 bool normalized = false);
265 /// return the approximation variance from each Approximation
266 /// (request forwarded to approxInterface)
267 const RealVector& approximation_variances(const Variables& vars);
268 /// return the approximation data from a particular Approximation
269 /// (request forwarded to approxInterface)
270 const Pecos::SurrogateData& approximation_data(size_t fn_index);
271
272 /// update component parallel mode for supporting parallelism in actualModel
273 void component_parallel_mode(short mode);
274
275 IntIntPair estimate_partition_bounds(int max_eval_concurrency);
276
277 /// set up actualModel for parallel operations
278 void derived_init_communicators(ParLevLIter pl_iter, int max_eval_concurrency,
279 bool recurse_flag = true);
280 /// set up actualModel for serial operations.
281 void derived_init_serial();
282 /// set active parallel configuration within actualModel
283 void derived_set_communicators(ParLevLIter pl_iter, int max_eval_concurrency,
284 bool recurse_flag = true);
285 /// deallocate communicator partitions for the DataFitSurrModel
286 /// (request forwarded to actualModel)
287 void derived_free_communicators(ParLevLIter pl_iter, int max_eval_concurrency,
288 bool recurse_flag = true);
289
290 /// Service actualModel job requests received from the master.
291 /// Completes when a termination message is received from stop_servers().
292 void serve_run(ParLevLIter pl_iter, int max_eval_concurrency);
293 /// Executed by the master to terminate actualModel server operations
294 /// when DataFitSurrModel iteration is complete.
295 void stop_servers();
296
297 /// update the Model's inactive view based on higher level (nested)
298 /// context and optionally recurse into actualModel
299 void inactive_view(short view, bool recurse_flag = true);
300
301 /// return the approxInterface identifier
302 const String& interface_id() const;
303 /// if recurse_flag, return the actualModel evaluation cache usage
304 bool evaluation_cache(bool recurse_flag = true) const;
305 /// if recurse_flag, return the actualModel restart file usage
306 bool restart_file(bool recurse_flag = true) const;
307
308 /// set the evaluation counter reference points for the DataFitSurrModel
309 /// (request forwarded to approxInterface and actualModel)
310 void set_evaluation_reference();
311 /// request fine-grained evaluation reporting within approxInterface
312 /// and actualModel
313 void fine_grained_evaluation_counters();
314 /// print the evaluation summary for the DataFitSurrModel
315 /// (request forwarded to approxInterface and actualModel)
316 void print_evaluation_summary(std::ostream& s, bool minimal_header = false,
317 bool relative_count = true) const;
318
319 /// set the warm start flag, including actualModel
320 void warm_start_flag(const bool flag);
321
322 ActiveSet default_interface_active_set();
323
324 //
325 //- Heading: Data members
326 //
327
328 /// whether to export the surrogate to file
329 const bool exportSurrogate;
330
331 /// whether to automatically refine the surrogate during the build phase
332 const bool autoRefine;
333 /// Maximum number of times to refine the surrogate
334 const int maxIterations;
335 /// Maximum number of evaluations while refining a surrogate
336 const int maxFuncEvals;
337 /// Convergence criterion, compared to CV score for specified metric
338 const Real convergenceTolerance;
339 /// Max number of iterations for which there is no average improvement
340 const int softConvergenceLimit;
341 /// Type of error metric to test for surrogate refinement convegence
342 const String refineCVMetric;
343 /// Number of cross validation folds for surrogate refinement
344 const int refineCVFolds;
345
346 private:
347
348 //
349 //- Heading: Convenience functions
350 //
351
352 /// optionally read surrogate data points from provided file
353 void import_points(unsigned short tabular_format, bool use_var_labels,
354 bool active_only);
355 /// initialize file stream for exporting surrogate evaluations
356 void initialize_export();
357 /// finalize file stream for exporting surrogate evaluations
358 void finalize_export();
359 /// initialize file stream for exporting surrogate evaluations
360 void export_point(int eval_id, const Variables& vars, const Response& resp);
361
362 /// Common code for processing of approximate response maps shared by
363 /// derived_synchronize() and derived_synchronize_nowait()
364 void derived_synchronize_approx(bool block,
365 IntResponseMap& approx_resp_map_rekey);
366
367 /// Updates fit arrays for local or multipoint approximations
368 void update_local_reference();
369 /// Builds a local or multipoint approximation using actualModel
370 void build_local_multipoint();
371 /// Builds a local or multipoint approximation using actualModel
372 void build_local_multipoint(const Variables& vars,
373 const IntResponsePair& response_pr);
374
375 /// Updates fit arrays for global approximations
376 void update_global_reference();
377 /// Builds a global approximation using daceIterator
378 void build_global();
379 /// Rebuilds a global approximation by generating new data using
380 /// daceIterator and appending to approxInterface
381 void rebuild_global();
382
383 /// Refine the built surrogate until convergence criteria are met
384 void refine_surrogate();
385
386 /// clear current data from approxInterface
387 void clear_approx_interface();
388 /// update anchor data in approxInterface
389 void update_approx_interface(const Variables& vars,
390 const IntResponsePair& response_pr);
391 /// build the approxInterface surrogate, passing variable bounds
392 void build_approx_interface();
393
394 /// update actualModel with data from constraints/labels/sets
395 void init_model(Model& model);
396 /// update actualModel with data from current variables/bounds
397 void update_model(Model& model);
398 /// update current variables/labels/bounds/targets with data from actualModel
399 void update_from_model(const Model& model);
400
401 /// test if inactive state is consistent
402 bool consistent(const Variables& vars) const;
403 /// test if active vars are within [l_bnds, u_bnds]
404 bool inside(const Variables& vars) const;
405 /// test for exact equality in values between active vars and sdv
406 bool active_vars_compare(const Variables& vars,
407 const Pecos::SurrogateDataVars& sdv) const;
408
409 //
410 //- Heading: Data members
411 //
412
413 /// manages construction and application of correction functions that
414 /// are applied to a surrogate model (DataFitSurr or HierarchSurr) in
415 /// order to reproduce high fidelity data.
416 DiscrepancyCorrection deltaCorr;
417
418 /// total points the user specified to construct the surrogate
419 int pointsTotal;
420 /// configuration for points management in build_global()
421 short pointsManagement;
422
423 /// type of point reuse for approximation builds: \c all, \c region
424 /// (default if points file), or \c none (default if no points file)
425 String pointReuse;
426 /// file name from \c import_build_points_file specification
427 String importPointsFile;
428
429 /// file name from \c export_approx_points_file specification
430 String exportPointsFile;
431 /// file export format for variables and approximate responses
432 unsigned short exportFormat;
433 /// output file stream for \c export_approx_points_file specification
434 std::ofstream exportFileStream;
435
436 /// file name from \c export_approx_variance_file specification
437 String exportVarianceFile;
438 /// file export format for variables and approximate response variance
439 unsigned short exportVarianceFormat;
440 /// output file stream for \c export_approx_variance_file specification
441 std::ofstream exportVarianceFileStream;
442
443 /// manages the building and subsequent evaluation of the approximations
444 /// (required for both global and local)
445 Interface approxInterface;
446
447 /// the truth model which provides evaluations for building the surrogate
448 /// (optional for global, required for local)
449 /** actualModel is unrestricted in type; arbitrary nestings are possible. */
450 Model actualModel;
451
452 /// selects parameter sets on which to evaluate actualModel in order
453 /// to generate the necessary data for building global approximations
454 /// (optional for global since restart data may also be used)
455 Iterator daceIterator;
456 };
457
458
~DataFitSurrModel()459 inline DataFitSurrModel::~DataFitSurrModel()
460 { if (!exportPointsFile.empty()) finalize_export(); }
461
462
463 inline void DataFitSurrModel::
nested_variable_mappings(const SizetArray & c_index1,const SizetArray & di_index1,const SizetArray & ds_index1,const SizetArray & dr_index1,const ShortArray & c_target2,const ShortArray & di_target2,const ShortArray & ds_target2,const ShortArray & dr_target2)464 nested_variable_mappings(const SizetArray& c_index1,
465 const SizetArray& di_index1,
466 const SizetArray& ds_index1,
467 const SizetArray& dr_index1,
468 const ShortArray& c_target2,
469 const ShortArray& di_target2,
470 const ShortArray& ds_target2,
471 const ShortArray& dr_target2)
472 {
473 // forward along to actualModel:
474 if (!actualModel.is_null())
475 actualModel.nested_variable_mappings(c_index1, di_index1, ds_index1,
476 dr_index1, c_target2, di_target2,
477 ds_target2, dr_target2);
478 }
479
480
nested_acv1_indices() const481 inline const SizetArray& DataFitSurrModel::nested_acv1_indices() const
482 { return actualModel.nested_acv1_indices(); }
483
484
nested_acv2_targets() const485 inline const ShortArray& DataFitSurrModel::nested_acv2_targets() const
486 { return actualModel.nested_acv2_targets(); }
487
488
query_distribution_parameter_derivatives() const489 inline short DataFitSurrModel::query_distribution_parameter_derivatives() const
490 {
491 return (actualModel.is_null()) ? NO_DERIVS :
492 actualModel.query_distribution_parameter_derivatives(); // forward along
493 }
494
495
qoi() const496 inline size_t DataFitSurrModel::qoi() const
497 {
498 switch (responseMode) {
499 // Response inflation from aggregation does not proliferate above
500 // this Model recursion level
501 case AGGREGATED_MODELS: return actualModel.qoi(); break;
502 default: return response_size(); break;
503 }
504 }
505
506
discrepancy_correction()507 inline DiscrepancyCorrection& DataFitSurrModel::discrepancy_correction()
508 { return deltaCorr; }
509
510
correction_type()511 inline short DataFitSurrModel::correction_type()
512 { return deltaCorr.correction_type(); }
513
514
correction_type(short corr_type)515 inline void DataFitSurrModel::correction_type(short corr_type)
516 { deltaCorr.correction_type(corr_type); }
517
518
total_points(int points)519 inline void DataFitSurrModel::total_points(int points)
520 { pointsTotal = points; if (points > 0) pointsManagement = TOTAL_POINTS; }
521
522
required_points()523 inline int DataFitSurrModel::required_points()
524 {
525 switch (pointsManagement) {
526 case TOTAL_POINTS: {
527 int min_points = approxInterface.minimum_points(true);
528 if (pointsTotal < min_points && outputLevel >= NORMAL_OUTPUT)
529 Cout << "\nDataFitSurrModel: Total points specified (" << pointsTotal
530 << ") is less than minimum required;\n "
531 << "increasing to " << min_points << std::endl;
532 return std::max(min_points, pointsTotal); break;
533 }
534 case RECOMMENDED_POINTS:
535 return approxInterface.recommended_points(true); break;
536 default: //case DEFAULT_POINTS: case MINIMUM_POINTS:
537 return approxInterface.minimum_points(true); break;
538 }
539 }
540
541
542 inline bool DataFitSurrModel::
active_vars_compare(const Variables & vars,const Pecos::SurrogateDataVars & sdv) const543 active_vars_compare(const Variables& vars,
544 const Pecos::SurrogateDataVars& sdv) const
545 {
546 // Similar to id_vars_exact_compare() in PRPMultiIndex.hpp
547
548 if (vars.is_null() || sdv.is_null())
549 return false;
550 // discrete strings not currently included in SurrogateDataVars
551 else if (vars.continuous_variables() != sdv.continuous_variables() ||
552 vars.discrete_int_variables() != sdv.discrete_int_variables() ||
553 vars.discrete_real_variables() != sdv.discrete_real_variables())
554 return false;
555
556 return true;
557 }
558
559
subordinate_iterator()560 inline Iterator& DataFitSurrModel::subordinate_iterator()
561 { return daceIterator; }
562
563
active_model_key(const UShortArray & key)564 inline void DataFitSurrModel::active_model_key(const UShortArray& key)
565 {
566 // assign activeKey and extract {surr,truth}ModelKey
567 SurrogateModel::active_model_key(key);
568
569 // recur both components: (actualModel could be hierarchical)
570 approxInterface.active_model_key(key);
571 actualModel.active_model_key(key);
572 }
573
574
clear_model_keys()575 inline void DataFitSurrModel::clear_model_keys()
576 { approxInterface.clear_model_keys(); }
577
578
surrogate_model()579 inline Model& DataFitSurrModel::surrogate_model()
580 {
581 // return by reference: OK to return letter instance
582 return *this;
583
584 // return by value: letter instance must be returned within an envelope for
585 // representation sharing/reference counting to work properly
586 //Model surr_model; // empty envelope
587 //surr_model.assign_rep(this); // populate letter, increment reference count
588 //return surr_model;
589 }
590
591
surrogate_model() const592 inline const Model& DataFitSurrModel::surrogate_model() const
593 { return *this; } // return of letter (see above)
594
595
truth_model()596 inline Model& DataFitSurrModel::truth_model()
597 { return actualModel; }
598
599
truth_model() const600 inline const Model& DataFitSurrModel::truth_model() const
601 { return actualModel; }
602
603
604 inline void DataFitSurrModel::
derived_subordinate_models(ModelList & ml,bool recurse_flag)605 derived_subordinate_models(ModelList& ml, bool recurse_flag)
606 {
607 if (!actualModel.is_null()) {
608 ml.push_back(actualModel);
609 if (recurse_flag)
610 actualModel.derived_subordinate_models(ml, recurse_flag);
611 }
612 }
613
614
resize_from_subordinate_model(size_t depth)615 inline void DataFitSurrModel::resize_from_subordinate_model(size_t depth)
616 {
617 if (!actualModel.is_null() && depth) {
618 // data flows from the bottom-up, so recurse first
619 if (depth == std::numeric_limits<size_t>::max())
620 actualModel.resize_from_subordinate_model(depth); // retain special value
621 else
622 actualModel.resize_from_subordinate_model(depth - 1);
623
624 // DataFitSurrModel consumes (newly) aggregated data sets through multiple
625 // SurrogateData instances --> don't resize locally.
626 //resize_response();
627 // It must therefore manage inflation of incoming ActiveSet instances
628
629 // daceIterator::activeSet muse be resized for consistency with actualModel
630 if (!daceIterator.is_null()) {
631 const ActiveSet& dace_set = daceIterator.active_set();
632 const ShortArray& dace_asv = dace_set.request_vector();
633 size_t num_am_resp = actualModel.response_size(),
634 num_dace_resp = dace_asv.size();
635 if (num_am_resp != num_dace_resp) {
636 ActiveSet new_set(dace_set); // deep copy
637 new_set.reshape(num_am_resp);
638 daceIterator.active_set(new_set);
639 }
640 }
641 }
642 }
643
644
update_from_subordinate_model(size_t depth)645 inline void DataFitSurrModel::update_from_subordinate_model(size_t depth)
646 {
647 if (!actualModel.is_null()) {
648 // data flows from the bottom-up, so recurse first
649 if (depth == std::numeric_limits<size_t>::max())
650 actualModel.update_from_subordinate_model(depth); // retain special value
651 else if (depth)
652 actualModel.update_from_subordinate_model(depth - 1);
653 // now pull the latest updates from actualModel
654 update_from_model(actualModel);
655 }
656 }
657
658
derived_interface()659 inline Interface& DataFitSurrModel::derived_interface()
660 { return approxInterface; }
661
662
663 inline void DataFitSurrModel::
primary_response_fn_weights(const RealVector & wts,bool recurse_flag)664 primary_response_fn_weights(const RealVector& wts, bool recurse_flag)
665 {
666 primaryRespFnWts = wts;
667 if (recurse_flag && !actualModel.is_null())
668 actualModel.primary_response_fn_weights(wts, recurse_flag);
669 }
670
671
surrogate_response_mode(short mode)672 inline void DataFitSurrModel::surrogate_response_mode(short mode)
673 {
674 responseMode = mode;
675
676 // Mode-specific logic:
677 // > Compared to HierarchSurrModel, don't need to be as strict in validating
678 // AUTO_CORRECTED_SURROGATE mode against corrType, since NO_CORRECTION is
679 // an admissible option in the case of global data fits. However,
680 // MODEL_DISCREPANCY still needs a discrepancy formulation (additive, etc.).
681 // > Management of multiple SurrogateData instances is complicated in the
682 // heterogeneous setting where level 0 uses a single instance and levels
683 // 1-L use two instances. In the future, could manage activation explicitly
684 // using functions shown below. For now, SurrogateData::{push,pop}() are
685 // hardened for inactive instances.
686 switch (mode) {
687 case MODEL_DISCREPANCY:
688 if (!corrType) {
689 Cerr << "Error: activation of mode MODEL_DISCREPANCY requires "
690 << "specification of a correction type." << std::endl;
691 abort_handler(MODEL_ERROR);
692 }
693 break;
694 case BYPASS_SURROGATE:
695 actualModel.surrogate_response_mode(mode); // recurse in this case
696 //approxInterface.deactivate_multilevel_approximation_data();
697 break;
698 case AGGREGATED_MODELS:
699 //approxInterface.activate_multilevel_approximation_data();
700 break;
701 }
702 }
703
704
705 //inline void DataFitSurrModel::link_multilevel_approximation_data()
706 //{ approxInterface.link_multilevel_approximation_data(); }
707
708
709 inline void DataFitSurrModel::
surrogate_function_indices(const IntSet & surr_fn_indices)710 surrogate_function_indices(const IntSet& surr_fn_indices)
711 {
712 surrogateFnIndices = surr_fn_indices;
713 approxInterface.approximation_function_indices(surr_fn_indices);
714 }
715
716
push_available()717 inline bool DataFitSurrModel::push_available()
718 { return approxInterface.push_available(); }
719
720
clear_inactive()721 inline void DataFitSurrModel::clear_inactive()
722 { approxInterface.clear_inactive(); }
723
724
advancement_available()725 inline bool DataFitSurrModel::advancement_available()
726 { return approxInterface.advancement_available(); }
727
728
formulation_updated() const729 inline bool DataFitSurrModel::formulation_updated() const
730 { return approxInterface.formulation_updated(); }
731
732
formulation_updated(bool update)733 inline void DataFitSurrModel::formulation_updated(bool update)
734 { approxInterface.formulation_updated(update); }
735
736
shared_approximation()737 inline SharedApproxData& DataFitSurrModel::shared_approximation()
738 { return approxInterface.shared_approximation(); }
739
740
approximations()741 inline std::vector<Approximation>& DataFitSurrModel::approximations()
742 { return approxInterface.approximations(); }
743
744
745 inline const RealVectorArray& DataFitSurrModel::
approximation_coefficients(bool normalized)746 approximation_coefficients(bool normalized)
747 { return approxInterface.approximation_coefficients(normalized); }
748
749
750 inline void DataFitSurrModel::
approximation_coefficients(const RealVectorArray & approx_coeffs,bool normalized)751 approximation_coefficients(const RealVectorArray& approx_coeffs,
752 bool normalized)
753 {
754 approxInterface.approximation_coefficients(approx_coeffs, normalized);
755
756 // Surrogate data is being imported. Update state to suppress automatic
757 // surrogate construction.
758 ++approxBuilds;
759 if (strbegins(surrogateType, "global_")) update_global_reference();
760 else update_local_reference();
761 }
762
763
764 inline void DataFitSurrModel::
rebuild_approximation(const IntResponsePair & response_pr)765 rebuild_approximation(const IntResponsePair& response_pr)
766 {
767 // decide which surrogates to rebuild based on resp_map content
768 BitArray rebuild_fns(numFns); // init to false
769 const ShortArray& asv = response_pr.second.active_set_request_vector();
770 for (size_t i=0; i<numFns; ++i)
771 if (asv[i])
772 rebuild_fns.set(i);
773 // rebuild the designated surrogates
774 approxInterface.rebuild_approximation(rebuild_fns);
775 ++approxBuilds;
776 }
777
778
779 inline void DataFitSurrModel::
rebuild_approximation(const IntResponseMap & all_resp)780 rebuild_approximation(const IntResponseMap& all_resp)
781 {
782 // decide which surrogates to rebuild based on resp_map content
783 BitArray rebuild_fns(numFns); // init to false
784 for (size_t i=0; i<numFns; ++i)
785 for (IntRespMCIter r_it=all_resp.begin(); r_it!=all_resp.end(); ++r_it)
786 if (r_it->second.active_set_request_vector()[i])
787 { rebuild_fns.set(i); break; }
788 // rebuild the designated surrogates
789 approxInterface.rebuild_approximation(rebuild_fns);
790 ++approxBuilds;
791 }
792
793
794 inline const RealVector& DataFitSurrModel::
approximation_variances(const Variables & vars)795 approximation_variances(const Variables& vars)
796 { return approxInterface.approximation_variances(vars); }
797
798
799 inline const Pecos::SurrogateData& DataFitSurrModel::
approximation_data(size_t fn_index)800 approximation_data(size_t fn_index)
801 { return approxInterface.approximation_data(fn_index); }
802
803
804 inline IntIntPair DataFitSurrModel::
estimate_partition_bounds(int max_eval_concurrency)805 estimate_partition_bounds(int max_eval_concurrency)
806 {
807 // support DB-based and on-the-fly instantiations for DataFitSurrModel
808 if (!daceIterator.is_null()) {
809 probDescDB.set_db_list_nodes(daceIterator.method_id());
810 return daceIterator.estimate_partition_bounds();
811 }
812 else if (!actualModel.is_null()) {
813 int am_max_conc = approxInterface.minimum_points(false)
814 * actualModel.derivative_concurrency(); // local/multipt
815 probDescDB.set_db_model_nodes(actualModel.model_id());
816 return actualModel.estimate_partition_bounds(am_max_conc);
817 }
818 else
819 return IntIntPair(1, 1);
820 }
821
822
derived_init_serial()823 inline void DataFitSurrModel::derived_init_serial()
824 {
825 //approxInterface.init_serial();
826
827 if (!actualModel.is_null())
828 actualModel.init_serial();
829 }
830
831
832 inline void DataFitSurrModel::
derived_set_communicators(ParLevLIter pl_iter,int max_eval_concurrency,bool recurse_flag)833 derived_set_communicators(ParLevLIter pl_iter, int max_eval_concurrency,
834 bool recurse_flag)
835 {
836 // allow recursion to progress - don't store/set/restore
837 //parallelLib.parallel_configuration_iterator(modelPCIter);
838 //approxInterface.set_communicators(messageLengths);
839 // asynchEvalFlag and evaluationCapacity updates not required for DFS
840 // (refer to {Recast,HierarchSurr}Model::derived_set_communicators())
841 //set_ie_asynchronous_mode(max_eval_concurrency);
842
843 miPLIndex = modelPCIter->mi_parallel_level_index(pl_iter);// run time setting
844
845 if (recurse_flag) {
846 if (!daceIterator.is_null())
847 daceIterator.set_communicators(pl_iter);
848 else if (!actualModel.is_null())
849 actualModel.init_communicators(pl_iter,
850 daceIterator.maximum_evaluation_concurrency()); // set in init_comms
851 }
852 }
853
854
855 inline void DataFitSurrModel::
derived_free_communicators(ParLevLIter pl_iter,int max_eval_concurrency,bool recurse_flag)856 derived_free_communicators(ParLevLIter pl_iter, int max_eval_concurrency,
857 bool recurse_flag)
858 {
859 // allow recursion to progress - don't store/set/restore
860 //parallelLib.parallel_configuration_iterator(modelPCIter);
861 //approxInterface.free_communicators();
862
863 if (recurse_flag) {
864 if (!daceIterator.is_null())
865 daceIterator.free_communicators(pl_iter);
866 else if (!actualModel.is_null())
867 actualModel.free_communicators(pl_iter,
868 daceIterator.maximum_evaluation_concurrency()); // set in init_comms
869 }
870 }
871
872
873 inline void DataFitSurrModel::
serve_run(ParLevLIter pl_iter,int max_eval_concurrency)874 serve_run(ParLevLIter pl_iter, int max_eval_concurrency)
875 {
876 // don't recurse, as actualModel.serve() will set actualModel comms
877 set_communicators(pl_iter, max_eval_concurrency, false);
878
879 if (!actualModel.is_null())
880 actualModel.serve_run(pl_iter,
881 daceIterator.maximum_evaluation_concurrency());
882 }
883
884
stop_servers()885 inline void DataFitSurrModel::stop_servers()
886 {
887 if (!actualModel.is_null())
888 actualModel.stop_servers();
889 }
890
891
inactive_view(short view,bool recurse_flag)892 inline void DataFitSurrModel::inactive_view(short view, bool recurse_flag)
893 {
894 currentVariables.inactive_view(view);
895 userDefinedConstraints.inactive_view(view);
896 if (recurse_flag && !actualModel.is_null())
897 actualModel.inactive_view(view, recurse_flag);
898 }
899
900
interface_id() const901 inline const String& DataFitSurrModel::interface_id() const
902 { return approxInterface.interface_id(); }
903
904
evaluation_cache(bool recurse_flag) const905 inline bool DataFitSurrModel::evaluation_cache(bool recurse_flag) const
906 {
907 return (recurse_flag && !actualModel.is_null()) ?
908 actualModel.evaluation_cache(recurse_flag) : false;
909 }
910
911
restart_file(bool recurse_flag) const912 inline bool DataFitSurrModel::restart_file(bool recurse_flag) const
913 {
914 return (recurse_flag && !actualModel.is_null()) ?
915 actualModel.restart_file(recurse_flag) : false;
916 }
917
918
set_evaluation_reference()919 inline void DataFitSurrModel::set_evaluation_reference()
920 {
921 approxInterface.set_evaluation_reference();
922
923 // don't recurse this, since the eval reference is for the top level iteration
924 //if (!actualModel.is_null())
925 // actualModel.set_evaluation_reference();
926
927 // may want to add this in time
928 //surrModelEvalRef = surrModelEvalCntr;
929 }
930
931
fine_grained_evaluation_counters()932 inline void DataFitSurrModel::fine_grained_evaluation_counters()
933 {
934 approxInterface.fine_grained_evaluation_counters(numFns);
935 if (!actualModel.is_null())
936 actualModel.fine_grained_evaluation_counters();
937 }
938
939
940 inline void DataFitSurrModel::
print_evaluation_summary(std::ostream & s,bool minimal_header,bool relative_count) const941 print_evaluation_summary(std::ostream& s, bool minimal_header,
942 bool relative_count) const
943 {
944 approxInterface.print_evaluation_summary(s, minimal_header, relative_count);
945 if (!actualModel.is_null()) {
946 if (daceIterator.is_null())
947 actualModel.print_evaluation_summary(s, minimal_header, relative_count);
948 else // daceIterator resets the eval reference -> don't use a relative count
949 actualModel.print_evaluation_summary(s, minimal_header, false);
950 }
951 }
952
953
warm_start_flag(const bool flag)954 inline void DataFitSurrModel::warm_start_flag(const bool flag)
955 {
956 warmStartFlag = flag;
957 actualModel.warm_start_flag(flag);
958 }
959
960 } // namespace Dakota
961
962 #endif
963