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:        DataModel
10 //- Description:
11 //-
12 //-
13 //- Owner:        Mike Eldred
14 //- Version: $Id: DataModel.hpp 6879 2010-07-30 01:05:11Z mseldre $
15 
16 #ifndef DATA_MODEL_H
17 #define DATA_MODEL_H
18 
19 #include "dakota_system_defs.hpp"
20 #include "dakota_data_types.hpp"
21 #include "MPIPackBuffer.hpp"
22 
23 namespace Dakota {
24 
25 /// define special values for pointsManagement
26 enum { DEFAULT_POINTS, MINIMUM_POINTS, RECOMMENDED_POINTS, TOTAL_POINTS };
27 
28 /// define special values for SurrogateModel::responseMode
29 enum { NO_SURROGATE=0,  UNCORRECTED_SURROGATE, AUTO_CORRECTED_SURROGATE,
30        BYPASS_SURROGATE, MODEL_DISCREPANCY, AGGREGATED_MODELS };
31 
32 /// define special values for approxCorrectionType
33 enum { NO_CORRECTION=0,  ADDITIVE_CORRECTION, MULTIPLICATIVE_CORRECTION,
34        COMBINED_CORRECTION };
35 
36 /// define types of random field approximations
37 enum { RF_KARHUNEN_LOEVE=0, RF_PCA_GP, RF_ICA };
38 
39 /// define types of analytic covariance functions
40 enum { NOCOVAR=0, EXP_L2, EXP_L1 };
41 
42 /// define special values for active subspace normalizations
43 enum { SUBSPACE_NORM_DEFAULT=0, SUBSPACE_NORM_MEAN_VALUE,
44        SUBSPACE_NORM_MEAN_GRAD, SUBSPACE_NORM_LOCAL_GRAD };
45 
46 /// define special values for componentParallelMode
47 /// (active model for parallel scheduling)
48 enum { NO_PARALLEL_MODE=0, SURROGATE_MODEL_MODE, TRUTH_MODEL_MODE,
49        SUB_MODEL_MODE, INTERFACE_MODE };
50 
51 /// define special values for distParamDerivs
52 enum { NO_DERIVS=0, ALL_DERIVS, MIXED_DERIVS };
53 
54 // define special values for regressionType in C3 FT (outside of Pecos).
55 // Note that C3 and Pecos are mutually exclusive: use of values from multiple
56 // enums should not conflict
57 enum { FT_LS, FT_RLS2 };//, FT_RLSD2, FT_RLSRKHS, FT_RLS1 };
58 // define special values for c3AdvanceType
59 enum { NO_C3_ADVANCEMENT=0, START_RANK_ADVANCEMENT, START_ORDER_ADVANCEMENT,
60        MAX_RANK_ADVANCEMENT, MAX_ORDER_ADVANCEMENT, MAX_RANK_ORDER_ADVANCEMENT};
61 
62 
63 /// Body class for model specification data.
64 
65 /** The DataModelRep class is used to contain the data from a model
66     keyword specification.  Default values are managed in the
67     DataModelRep constructor.  Data is public to avoid maintaining
68     set/get functions, but is still encapsulated within ProblemDescDB
69     since ProblemDescDB::dataModelList is private. */
70 
71 class DataModelRep
72 {
73   //
74   //- Heading: Friends
75   //
76 
77   /// the handle class can access attributes of the body class directly
78   friend class DataModel;
79 
80 public:
81 
82   ~DataModelRep(); ///< destructor
83 
84   //
85   //- Heading: Data
86   //
87 
88   /// string identifier for the model specification data set (from
89   /// the \c id_model specification in \ref ModelIndControl)
90   String idModel;
91   /// model type selection: single, surrogate, or nested (from the model type
92   /// specification in \ref ModelIndControl)
93   String modelType;
94   /// string pointer to the variables specification to be used by this model
95   /// (from the \c variables_pointer specification in \ref ModelIndControl)
96   String variablesPointer;
97   /// string pointer to the interface specification to be used by this model
98   /// (from the \c interface_pointer specification in \ref ModelSingle and
99   /// the \c optional_interface_pointer specification in \ref ModelNested)
100   String interfacePointer;
101   /// string pointer to the responses specification to be used by this model
102   /// (from the \c responses_pointer specification in \ref ModelIndControl)
103   String responsesPointer;
104   /// whether this model and its children will add hierarchy-based
105   /// tags to eval ids
106   bool hierarchicalTags;
107   /// pointer to a sub-iterator used for global approximations (from the
108   /// \c dace_method_pointer specification in \ref ModelSurrG) or by
109   /// nested models (from the \c sub_method_pointer specification in
110   /// \ref ModelNested)
111   String subMethodPointer;
112 
113   // single/simulation models
114 
115   /// (state) variable identifier that defines a set or range of solution
116   /// level controls (space/time discretization levels, iterative convergence
117   /// tolerances, etc.) for defining a secondary hierarchy of fidelity within
118   /// the scope of a single model form (from \c solution_level_control
119   /// specification; see also \c ordered_model_fidelities)
120   String solutionLevelControl;
121   /// array of relative simulation costs corresponding to each of the
122   /// solution levels (from \c solution_level_cost specification; see also
123   /// \c solution_level_control); a scalar input is interpreted as a constant
124   /// cost multiplier to be applied recursively
125   RealVector solutionLevelCost;
126 
127   // surrogate models
128 
129   /// array specifying the response function set that is approximated
130   IntSet surrogateFnIndices;
131   /// the selected surrogate type: local_taylor, multipoint_tana,
132   /// global_(neural_network,mars,orthogonal_polynomial,gaussian,
133   /// polynomial,kriging), or hierarchical
134   String surrogateType;
135   /// pointer to the model specification for constructing the truth model
136   /// used in constructing surrogates (from the \c actual_model_pointer
137   /// specification in \ref ModelSurrL and \ref ModelSurrMP)
138   String actualModelPointer;
139   /// an ordered list of model pointers (low to high) corresponding to a
140   /// hierarchy of modeling fidelity (from the \c ordered_model_fidelities
141   /// specification in \ref ModelSurrH)
142   StringArray orderedModelPointers;
143 
144   // controls for number of points with which to build the model
145 
146   /// user-specified lower bound on total points with which to build the model
147   /// (if reuse_points < pointsTotal, new samples will make up the difference)
148   int pointsTotal;
149   /// points management configuration for DataFitSurrModel:
150   /// DEFAULT_POINTS, MINIMUM_POINTS, or RECOMMENDED_POINTS
151   short pointsManagement;
152 
153   /// sample reuse selection for building global approximations: none, all,
154   /// region, or file (from the \c reuse_samples specification in
155   /// \ref ModelSurrG)
156   String approxPointReuse;
157 
158   /// the file name from the \c import_build_points_file specification in
159   /// \ref ModelSurrG
160   String importBuildPtsFile;
161   /// tabular format for the build point import file
162   unsigned short importBuildFormat;
163   /// whether to parse/validate variable labels from header
164   bool importUseVariableLabels;
165   /// whether to import active variables only
166   bool importBuildActive;
167 
168   // the file name from the \c import_approx_points_file specification in
169   // \ref ModelSurrG
170   //String importApproxPtsFile;
171   // tabular format for the approx point import file
172   //unsigned short importApproxFormat;
173   // whether to import active variables only
174   //bool importApproxActive;
175 
176   /// the file name from the \c export_approx_points_file specification in
177   /// \ref ModelSurrG
178   String exportApproxPtsFile;
179   /// tabular format for the approx point export file
180   unsigned short exportApproxFormat;
181 
182   /// filename for surrogate variance evaluation export
183   String exportApproxVarianceFile;
184   /// tabular format for the approx variance export file
185   unsigned short exportApproxVarianceFormat;
186 
187   /// Option to turn on surrogate model export (export_model)
188   bool exportSurrogate;
189 
190   /// the filename prefix for export_model
191   String modelExportPrefix;
192 
193   /// Format selection for export_model
194   unsigned short modelExportFormat;
195 
196   /// correction type for global and hierarchical approximations:
197   /// NO_CORRECTION, ADDITIVE_CORRECTION, MULTIPLICATIVE_CORRECTION,
198   /// or COMBINED_CORRECTION (from the \c correction specification in
199   /// \ref ModelSurrG and \ref ModelSurrH)
200   short approxCorrectionType;
201   /// correction order for global and hierarchical approximations: 0,
202   /// 1, or 2 (from the \c correction specification in \ref ModelSurrG
203   /// and \ref ModelSurrH)
204   short approxCorrectionOrder;
205   /// flags the use of derivatives in building global approximations
206   /// (from the \c use_derivatives specification in \ref ModelSurrG)
207   bool modelUseDerivsFlag;
208   /// scalar integer indicating the order of the polynomial approximation
209   /// (1=linear, 2=quadratic, 3=cubic; from the \c polynomial specification
210   /// in \ref ModelSurrG)
211   short polynomialOrder;
212   /// vector of correlations used in building a kriging approximation
213   /// (from the \c correlations specification in \ref ModelSurrG)
214   RealVector krigingCorrelations;
215   /// optimization method to use in finding optimal correlation parameters:
216   /// none, sampling, local, global
217   String krigingOptMethod;
218   /// maximum number of trials in optimization of kriging correlations
219   short krigingMaxTrials;
220   /// upper bound on kriging correlation vector
221   RealVector krigingMaxCorrelations;
222   /// lower bound on kriging correlation vector
223   RealVector krigingMinCorrelations;
224   /// nugget value for kriging
225   Real krigingNugget;
226   /// option to have Kriging find the best nugget value to use
227   short krigingFindNugget;
228   /// weight function for moving least squares approximation
229   short mlsWeightFunction;
230   /// bases for radial basis function approximation
231   short rbfBases;
232   /// maximum number of points for radial basis function approximation
233   short rbfMaxPts;
234   /// maximum number of subsets for radial basis function approximation
235   short rbfMaxSubsets;
236   /// minimum partition for radial basis function approximation
237   short rbfMinPartition;
238   /// maximum number of bases for MARS approximation
239   short marsMaxBases;
240   /// interpolation type for MARS approximation
241   String marsInterpolation;
242   /// random weight for artificial neural network approximation
243   short annRandomWeight;
244   /// number of nodes for artificial neural network approximation
245   short annNodes;
246   /// range for artificial neural network approximation
247   Real annRange;
248   /// number of restarts for gradient-based optimization in GP
249   int numRestarts;
250 
251   /// whether domain decomposition is enabled
252   bool domainDecomp;
253   /// type of local cell of domain decomp
254   String decompCellType;
255   /// number of support layers for each local basis function
256   int decompSupportLayers;
257   /// whether discontinuity detection is enabled
258   bool decompDiscontDetect;
259   /// function value (jump) threshold for discontinuity detection in
260   /// domain decomp
261   Real discontJumpThresh;
262   /// gradient threshold for discontinuity detection in domain decomp
263   Real discontGradThresh;
264 
265   /// scalar integer indicating the order of the Gaussian process mean
266   /// (0= constant, 1=linear, 2=quadratic, 3=cubic); from the
267   /// \c gaussian_process specification  in \ref ModelSurrG)
268   String trendOrder;
269   /// flag indicating the use of point selection in the Gaussian process
270   bool pointSelection;
271   /// List of diagnostic metrics the user requests to assess the
272   /// goodness of fit for a surrogate model.
273   StringArray diagMetrics;
274   /// flag indicating the use of cross validation on the metrics specified
275   bool crossValidateFlag;
276   /// number of folds to perform in cross validation
277   int numFolds;
278   /// percentage of data to withhold for cross validation process
279   Real percentFold;
280   /// flag indicating the use of PRESS on the metrics specified
281   bool pressFlag;
282 
283   /// the file name from the \c challenge_points_file specification in
284   /// \ref ModelSurrG
285   String importChallengePtsFile;
286   /// tabular format of the challenge data file
287   unsigned short importChallengeFormat;
288   /// whether to parse/validate variable labels from header
289   bool importChalUseVariableLabels;
290   /// whether to import active variables only
291   bool importChallengeActive;
292 
293   /// file containing advanced surrogate option overrides
294   String advancedOptionsFilename;
295 
296   // nested models
297 
298   /// string pointer to the responses specification used by the optional
299   /// interface in nested models (from the \c
300   /// optional_interface_responses_pointer specification in \ref ModelNested)
301   String optionalInterfRespPointer;
302   /// the primary variable mappings used in nested models for identifying
303   /// the lower level variable targets for inserting top level variable
304   /// values (from the \c primary_variable_mapping specification in
305   /// \ref ModelNested)
306   StringArray primaryVarMaps;
307   /// the secondary variable mappings used in nested models for identifying
308   /// the (distribution) parameter targets within the lower level variables
309   /// for inserting top level variable values (from the
310   /// \c secondary_variable_mapping specification in \ref ModelNested)
311   StringArray secondaryVarMaps;
312   /// the primary response mapping matrix used in nested models for
313   /// weighting contributions from the sub-iterator responses in the top
314   /// level (objective) functions (from the \c primary_response_mapping
315   /// specification in \ref ModelNested)
316   RealVector primaryRespCoeffs;
317   /// the secondary response mapping matrix used in nested models for
318   /// weighting contributions from the sub-iterator responses in the top
319   /// level (constraint) functions (from the \c secondary_response_mapping
320   /// specification in \ref ModelNested)
321   RealVector secondaryRespCoeffs;
322   /// whether an identity response map is requested in lieu of explicit maps
323   bool identityRespMap;
324   /// number of servers for concurrent sub-iterator parallelism
325   int subMethodServers;
326   /// number of processors for each concurrent sub-iterator partition
327   int subMethodProcs;
328   /// scheduling approach for concurrent sub-iterator parallelism:
329   /// {DEFAULT,MASTER,PEER}_SCHEDULING
330   short subMethodScheduling;
331 
332   // subspace models
333 
334   /// initial samples to build the subspace model
335   int initialSamples;
336 
337   /// sampling method for building the subspace model
338   unsigned short subspaceSampleType;
339 
340   /// refinement samples to add in each batch
341   IntVector refineSamples;
342   /// maximum number of subspace build iterations
343   int maxIterations;
344   /// convergence tolerance on build process
345   Real convergenceTolerance;
346 
347   /// Flag to use Bing Li method to identify active subspace dimension
348   bool subspaceIdBingLi;
349 
350   /// Flag to use Constantine method to identify active subspace dimension
351   bool subspaceIdConstantine;
352 
353   /// Flag to use eigenvalue energy method to identify active subspace dimension
354   bool subspaceIdEnergy;
355 
356   /// Flag to build surrogate over active subspace
357   bool subspaceBuildSurrogate;
358 
359   /// Size of subspace
360   int subspaceDimension;
361 
362   /// Normalization to use when forming a subspace with multiple response
363   /// functions
364   unsigned short subspaceNormalization;
365 
366   /// Number of bootstrap samples for subspace identification
367   int numReplicates;
368 
369   /// Flag to use cross validation to identify active subspace dimension
370   bool subspaceIdCV;
371 
372   /// relative tolerance used by cross validation subspace dimension id method
373   Real relTolerance;
374 
375   /// decrease tolerance used by cross validation subspace dimension id method
376   Real decreaseTolerance;
377 
378   /// maximum rank considered by cross validation subspace dimension id method
379   int subspaceCVMaxRank;
380 
381   /// flag to use incremental dimension estimation in the cross validation metric
382   bool subspaceCVIncremental;
383 
384   /// Contains which cutoff method to use in the cross validation metric
385   unsigned short subspaceIdCVMethod;
386 
387   // Function-Train Options
388 
389   /// type of (regularized) regression: FT_LS or FT_RLS2
390   short regressionType;
391   /// penalty parameter for regularized regression (FT_RLS2)
392   Real regressionL2Penalty;
393   /// max iterations for optimization solver used in FT regression
394   int maxSolverIterations;
395   /// maximum number of cross iterations
396   int maxCrossIterations;
397   /// optimization tolerance for FT regression
398   Real solverTol;
399   /// Rounding tolerance for FT regression
400   Real solverRoundingTol;
401   /// arithmetic (rounding) tolerance for FT sums and products
402   Real statsRoundingTol;
403   /// sub-sample a tensor grid for generating regression data
404   bool tensorGridFlag;
405   /// starting polynomial order
406   unsigned short startOrder;
407   /// polynomial order increment when adapting
408   unsigned short kickOrder;
409   /// maximum order of basis polynomials
410   unsigned short maxOrder;
411   /// whether or not to adapt order by cross validation
412   bool adaptOrder;
413   /// starting rank
414   size_t startRank;
415   /// rank increase increment
416   size_t kickRank;
417   /// maximum rank
418   size_t maxRank;
419   /// whether or not to adapt rank
420   bool adaptRank;
421   /// quantity to increment (start rank, start order, max rank, max order,
422   /// max rank + max order) for FT (uniform) p-refinement
423   short c3AdvanceType;
424   // refinement type for stochastic expansions: P_REFINEMENT, H_REFINEMENT
425   //short refinementType;
426   // refinement control for stochastic expansions: UNIFORM, DIMENSION_ADAPTIVE
427   //short refinementControl;
428 
429   /// number of data points used in FT construction by regression
430   size_t collocationPoints;
431   /// ratio of number of points to nuber of unknowns
432   Real collocationRatio;
433 
434   /// whether automatic surrogate refinement is enabled
435   bool autoRefine;
436   /// maximum evals in refinement
437   int maxFunctionEvals;
438   /// metric to use in cross-validation guided refinement
439   String refineCVMetric;
440   /// max number of iterations in refinement without improvement
441   int softConvergenceLimit;
442   /// number of cross-validation folds in guided refinement
443   int refineCVFolds;
444 
445   /// sparse grid level for low-order PCE used to compute rotation
446   /// matrix within adapted basis approach to dimension reduction
447   unsigned short adaptedBasisSparseGridLev;
448   /// expansion order for low-order PCE used to compute rotation
449   /// matrix within adapted basis approach to dimension reduction
450   unsigned short adaptedBasisExpOrder;
451   /// collocation ratio for low-order PCE used to compute rotation
452   /// matrix within adapted basis approach to dimension reduction
453   Real adaptedBasisCollocRatio;
454 
455   // random field models
456 
457   /// Contains which type of random field model
458   unsigned short randomFieldIdForm;
459 
460   /// Contains which type of analytic covariance function
461   unsigned short analyticCovIdForm;
462 
463   /// truncation tolerance on build process: percent variance explained
464   Real truncationTolerance;
465 
466   /// pointer to the model through which to propagate the random field
467   String propagationModelPointer;
468   /// File from which to build the random field
469   String rfDataFileName;
470 
471 private:
472 
473   //
474   //- Heading: Constructors, destructor, operators
475   //
476 
477   DataModelRep();  ///< constructor
478 
479   //
480   //- Heading: Member methods
481   //
482 
483   /// write a DataModelRep object to an std::ostream
484   void write(std::ostream& s) const;
485 
486   /// read a DataModelRep object from a packed MPI buffer
487   void read(MPIUnpackBuffer& s);
488   /// write a DataModelRep object to a packed MPI buffer
489   void write(MPIPackBuffer& s) const;
490 
491   //
492   //- Heading: Private convenience functions
493   //
494 
495 };
496 
497 
~DataModelRep()498 inline DataModelRep::~DataModelRep() { }
499 
500 
501 /// Handle class for model specification data.
502 
503 /** The DataModel class is used to provide a memory management handle
504     for the data in DataModelRep.  It is populated by
505     IDRProblemDescDB::model_kwhandler() and is queried by the
506     ProblemDescDB::get_<datatype>() functions.  A list of DataModel
507     objects is maintained in ProblemDescDB::dataModelList, one for
508     each model specification in an input file. */
509 
510 class DataModel
511 {
512   //
513   //- Heading: Friends
514   //
515 
516   // the problem description database
517   friend class ProblemDescDB;
518   // the NIDR derived problem description database
519   friend class NIDRProblemDescDB;
520 
521 public:
522 
523   /// compares the idModel attribute of DataModel objects
id_compare(const DataModel & dm,const std::string & id)524   static bool id_compare(const DataModel& dm, const std::string& id)
525   { return id == dm.dataModelRep->idModel; }
526 
527   //
528   //- Heading: Constructors, destructor, operators
529   //
530 
531   DataModel();                                ///< constructor
532   DataModel(const DataModel&);               ///< copy constructor
533   ~DataModel();                               ///< destructor
534 
535   DataModel& operator=(const DataModel&); ///< assignment operator
536 
537   //
538   //- Heading: Member methods
539   //
540 
541   /// write a DataModel object to an std::ostream
542   void write(std::ostream& s) const;
543 
544   /// read a DataModel object from a packed MPI buffer
545   void read(MPIUnpackBuffer& s);
546   /// write a DataModel object to a packed MPI buffer
547   void write(MPIPackBuffer& s) const;
548 
549   /// return dataModelRep
550   std::shared_ptr<DataModelRep> data_rep();
551 
552 private:
553 
554   //
555   //- Heading: Data
556   //
557 
558   /// pointer to the body (handle-body idiom)
559   std::shared_ptr<DataModelRep> dataModelRep;
560 };
561 
562 
data_rep()563 inline std::shared_ptr<DataModelRep> DataModel::data_rep()
564 {return dataModelRep; }
565 
566 
567 /// MPIPackBuffer insertion operator for DataModel
operator <<(MPIPackBuffer & s,const DataModel & data)568 inline MPIPackBuffer& operator<<(MPIPackBuffer& s, const DataModel& data)
569 { data.write(s); return s; }
570 
571 
572 /// MPIUnpackBuffer extraction operator for DataModel
operator >>(MPIUnpackBuffer & s,DataModel & data)573 inline MPIUnpackBuffer& operator>>(MPIUnpackBuffer& s, DataModel& data)
574 { data.read(s); return s; }
575 
576 
577 /// std::ostream insertion operator for DataModel
operator <<(std::ostream & s,const DataModel & data)578 inline std::ostream& operator<<(std::ostream& s, const DataModel& data)
579 { data.write(s); return s; }
580 
write(std::ostream & s) const581 inline void DataModel::write(std::ostream& s) const
582 { dataModelRep->write(s); }
583 
584 
read(MPIUnpackBuffer & s)585 inline void DataModel::read(MPIUnpackBuffer& s)
586 { dataModelRep->read(s); }
587 
588 
write(MPIPackBuffer & s) const589 inline void DataModel::write(MPIPackBuffer& s) const
590 { dataModelRep->write(s); }
591 
592 } // namespace Dakota
593 
594 #endif
595