[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
base_parameterization.cpp
1 #include "base_parameterization.hpp"
2 
3 namespace PHiLiP {
4 
5 template<int dim>
7  std::shared_ptr<HighOrderGrid<dim,double>> _high_order_grid)
8  : high_order_grid(_high_order_grid)
9  , mpi_communicator(MPI_COMM_WORLD)
10  , pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(mpi_communicator)==0)
11 {
12  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
13  MPI_Comm_size(MPI_COMM_WORLD, &n_mpi);
14 }
15 
16 template<int dim>
17 void BaseParameterization<dim> :: output_design_variables(unsigned int /*iteration_no*/) const
18 {
19  // Outputs nothing by default. Overriden in derived classes.
20 }
21 
22 template<int dim>
24  const VectorType &current_design_var,
25  const VectorType &updated_design_var) const
26 {
27  VectorType diff = current_design_var;
28  diff -= updated_design_var;
29 
30  bool is_design_variable_updated = true;
31  if(diff.l2_norm() == 0.0) {is_design_variable_updated = false;}
32 
33  return is_design_variable_updated;
34 }
35 
37 } // PHiLiP namespace
int n_mpi
Total no. of processors.
int mpi_rank
Processor# of current processor.
BaseParameterization(std::shared_ptr< HighOrderGrid< dim, double >> _high_order_grid)
Constructor.
bool has_design_variable_been_updated(const VectorType &current_design_var, const VectorType &updated_design_var) const
Checks if the design variable has changed.
Files for the baseline physics.
Definition: ADTypes.hpp:10
Abstract class for design parameterization. Objective function and the constraints take this class&#39;s ...
dealii::LinearAlgebra::distributed::Vector< double > VectorType
Alias for dealii&#39;s parallel distributed vector.
virtual void output_design_variables(const unsigned int) const
Outputs design variables. Doesn&#39;t output anything if not overridden.