[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
halton_sampling.cpp
1 #include "halton_sampling.h"
2 #include <iostream>
3 #include <filesystem>
4 #include "functional/functional.h"
5 #include <deal.II/lac/trilinos_sparse_matrix.h>
6 #include <deal.II/lac/vector_operation.h>
7 #include "reduced_order_solution.h"
8 #include "flow_solver/flow_solver.h"
9 #include "flow_solver/flow_solver_factory.h"
10 #include <cmath>
11 #include "rbf_interpolation.h"
12 #include "ROL_Algorithm.hpp"
13 #include "ROL_LineSearchStep.hpp"
14 #include "ROL_StatusTest.hpp"
15 #include "ROL_Stream.hpp"
16 #include "ROL_Bounds.hpp"
17 #include "halton.h"
18 #include "min_max_scaler.h"
19 #include <deal.II/base/timer.h>
20 
21 namespace PHiLiP {
22 
23 template<int dim, int nstate>
25  const dealii::ParameterHandler &parameter_handler_input)
26  : AdaptiveSamplingBase<dim, nstate>(parameters_input, parameter_handler_input)
27 {}
28 
29 template <int dim, int nstate>
31 {
32  this->pcout << "Starting Halton sampling process" << std::endl;
33  auto stream = this->pcout;
34  dealii::TimerOutput timer(stream,dealii::TimerOutput::summary,dealii::TimerOutput::wall_times);
35  timer.enter_subsection ("Solve FOMs and Assemble POD");
36  this->placeInitialSnapshots();
37  this->current_pod->computeBasis();
38 
39  timer.leave_subsection();
40 
41  this->outputIterationData("final");
42 
43  return 0;
44 }
45 
46 #if PHILIP_DIM==1
48 #endif
49 
50 #if PHILIP_DIM!=1
52 #endif
53 
54 }
dealii::ConditionalOStream pcout
ConditionalOStream.
Files for the baseline physics.
Definition: ADTypes.hpp:10
Halton sampling.
Main parameter class that contains the various other sub-parameter classes.
HaltonSampling(const PHiLiP::Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Constructor.
std::shared_ptr< ProperOrthogonalDecomposition::OnlinePOD< dim > > current_pod
Most up to date POD basis.
virtual void outputIterationData(std::string iteration) const
Output for each iteration.
void placeInitialSnapshots() const
Placement of initial snapshots.
int run_sampling() const override
Run Sampling Procedure.