17 #ifndef _FUEL_CELL_APPLICATION_CORE_OPTIMIZATION_BLOCK_MATRIX_APPLICATION_H_ 
   18 #define _FUEL_CELL_APPLICATION_CORE_OPTIMIZATION_BLOCK_MATRIX_APPLICATION_H_ 
   21 #include <deal.II/base/convergence_table.h> 
   22 #include <deal.II/grid/persistent_tria.h> 
   23 #include <deal.II/dofs/dof_renumbering.h> 
   32 namespace ApplicationCore
 
   75                                                bool triangulation_only);
 
  104             virtual void initialize(ParameterHandler& param);
 
  112             virtual void responses (std::vector<double>& f,
 
  161             virtual void dresponses_dl (std::vector<std::vector<double> >& df_dl,
 
  193             virtual void dresponses_du(std::vector<FuelCell::ApplicationCore::FEVector>& dst,
 
  200             virtual void cell_dresponses_du(std::vector<FuelCell::ApplicationCore::FEVector>& df_du,
 
  202                                             std::vector<std::vector<double> >& src);
 
  227             virtual void dresidual_dlambda(std::vector<FuelCell::ApplicationCore::FEVector>& dst,
 
  235                                                 std::vector<std::vector<double> >& src);
 
  248             void solve_direct (std::vector<std::vector<double> >& df_dl,
 
  262             void solve_adjoint (std::vector<std::vector<double> >& df_dl,
 
  311                 ParameterHandler prm;
 
  314                 file.open(
"Default_parameter_file.prm");
 
  315                 prm.print_parameters (file, ParameterHandler::Text);
 
  353                                                                   const unsigned int&              refinement_cycle,
 
  354                                                                   std::vector< ConvergenceTable >& convergence_tables)
 const 
  356                 const std::type_info& info = 
typeid(*this);
 
  357                 FcstUtilities::log << 
"Pure function " << __FUNCTION__ << 
" called in Class " << info.name() << std::endl;
 
  368                 for (
unsigned int i=0; i<
n_resp; ++i)
 
  369                     for (
unsigned int j=0; j<
n_dvar; ++j)
 
  381                 for (
unsigned int r=0; r<
n_resp; ++r)
 
  384                     unsigned int n_blocks = df_du[r].n_blocks();
 
  385                     unsigned int size = df_du[r].size();
 
  386                     unsigned int comp_per_block = size/n_blocks;
 
  387                     for (
unsigned int i=0; i<n_blocks; ++i)
 
  388                         for (
unsigned int j=0; j<comp_per_block; ++j)
 
  390                             FcstUtilities::log<<
"Element "<<j<<
" in block "<<i<<
" in STEP (Du) is "<<df_du[r].block(i)(j)<<std::endl;
 
  402             const std::string 
extend_filename(
const std::string& , 
const int precision = 3) 
const;
 
  409             typedef void (*
CELL_Dvalues) (std::vector<FuelCell::ApplicationCore::FEVector>& ,
 
  411                                           std::vector<std::vector<double> >& );
 
  422             void dfunction (std::vector<FuelCell::ApplicationCore::FEVector>& dst,
 
void solve_adjoint(std::vector< std::vector< double > > &df_dl, const FuelCell::ApplicationCore::FEVector &sol)
Solver in order to obtain the analytical sensitivities for the given objective and constraints using ...
 
virtual void dresidual_dlambda(std::vector< FuelCell::ApplicationCore::FEVector > &dst, const FuelCell::ApplicationCore::FEVectors &src)
Loop over all cells and assemble the vector  that is used to solve the sensitivity equations by using...
 
void dfunction(std::vector< FuelCell::ApplicationCore::FEVector > &dst, const FuelCell::ApplicationCore::FEVectors &src, bool dfunctional_du, bool dresidual_dlambda)
This is an auxiliary function called by dresidual_dlambda and dfunctional_du. 
 
virtual void dresponses_du(std::vector< FuelCell::ApplicationCore::FEVector > &dst, const FuelCell::ApplicationCore::FEVectors &src)
Loop over all cells and assemble the vector  that is used to solve the sensitivity equations by using...
 
void solve_direct(std::vector< std::vector< double > > &df_dl, const FuelCell::ApplicationCore::FEVectors &sol)
Solver in order to obtain the analytical sensitivities for the given objective and constraints using ...
 
unsigned int n_obj
Number of objective functions. 
Definition: optimization_block_matrix_application.h:440
 
void set_output_variables(std::vector< std::string > &dakota_name_responses)
Definition: optimization_block_matrix_application.h:334
 
virtual void dresponses_dl(std::vector< std::vector< double > > &df_dl, const FuelCell::ApplicationCore::FEVectors &src)
Post-processing. 
 
Application handling matrices and assembling linear systems of equations. 
Definition: block_matrix_application.h:74
 
virtual void global_dresponses_du(std::vector< FuelCell::ApplicationCore::FEVector > &df_du, const FuelCell::ApplicationCore::FEVector &src)
This class is used to evaluate the sensitivities of all responses that do not require looping over ce...
 
void(* CELL_Dvalues)(std::vector< FuelCell::ApplicationCore::FEVector > &, const typename DoFApplication< dim >::CellInfo &, std::vector< std::vector< double > > &)
Definition of a pointer to a function that is called in a loop over cells to compute the derivatives ...
Definition: optimization_block_matrix_application.h:409
 
virtual void cell_dresponses_du(std::vector< FuelCell::ApplicationCore::FEVector > &df_du, const typename DoFApplication< dim >::CellInfo &info, std::vector< std::vector< double > > &src)
Integration of local bilinear form. 
 
virtual void global_dresponses_dl(std::vector< std::vector< double > > &df_dl, const FuelCell::ApplicationCore::FEVector &src)
This class is used to evaluate the sensitivities of all responses that do not require looping over ce...
 
std::vector< std::string > name_responses
Member that stores the name of the responses, i.e. 
Definition: optimization_block_matrix_application.h:453
 
OptimizationBlockMatrixApplication(FuelCell::ApplicationCore::DoFApplication< dim > &, bool triangulation_only)
Constructor for an object, owning its own mesh and dof handler. 
 
virtual void initialize(ParameterHandler ¶m)
Initialize application. 
 
std::vector< unsigned int > user_input_bdry
Boundary id on which the boundary response is computed. 
Definition: optimization_block_matrix_application.h:480
 
This class is created for the objects handed to the mesh loops. 
Definition: mesh_loop_info_objects.h:544
 
std::vector< std::string > all_response_names
Variable that holds a list of response names. 
Definition: optimization_block_matrix_application.h:475
 
std::vector< std::string > name_design_var
Member that stores the name of the design variables. 
Definition: optimization_block_matrix_application.h:448
 
unsigned int get_n_resp() const 
Member function that returns the number of responses. 
 
void set_optimization_parameters(unsigned int &n_dvar, unsigned int &n_resp, std::vector< std::string > &name_design_var, std::vector< std::string > &name_responses)
This routine assigns the value of n_dvar, n_resp, name_design_var and name_responses to the internal ...
 
virtual void print_responses(std::vector< double > &resp)
This function is used to print the responses to screen (FcstUtilities::log) 
 
virtual void cell_dresidual_dlambda(std::vector< FuelCell::ApplicationCore::FEVector > &cell_vector, const typename DoFApplication< dim >::CellInfo &cell, std::vector< std::vector< double > > &src)
Integration of local bilinear form. 
 
virtual void compute_L1_L2_error_and_convergence_rate(const FuelCell::ApplicationCore::FEVector &solution, const unsigned int &refinement_cycle, std::vector< ConvergenceTable > &convergence_tables) const 
If the exact or analytical solution is available, then this function computes L1 and L2 norms of the ...
Definition: optimization_block_matrix_application.h:352
 
bool output_coarse_solution
Decision variable for whether the solution is to be output for transfer. 
Definition: optimization_block_matrix_application.h:462
 
bool get_bool_transfer_solution()
Function to see if we are transferring a solution on a refined grid to the initial coarse grid...
Definition: optimization_block_matrix_application.h:320
 
FCSTLogStream log
Object used to output data to file and, if file attached recorded to a file as well. 
 
virtual void responses(std::vector< double > &f, const FuelCell::ApplicationCore::FEVectors &vectors)
Post-processing. 
 
void print_dresponses_du(std::vector< FuelCell::ApplicationCore::FEVector > df_du)
Auxiliary routine to print the values of df_du This routine should be called once df_du is assembled...
Definition: optimization_block_matrix_application.h:377
 
void print_default_parameter_file()
Print the default parameter handler file. 
Definition: optimization_block_matrix_application.h:309
 
bool boundary_responses
true, if boundary responses are supposed to be computed. 
Definition: optimization_block_matrix_application.h:470
 
std::vector< std::string > name_output_var
Member that stores the name of the output variables, These names will be written to name_responses if...
Definition: optimization_block_matrix_application.h:459
 
virtual void check_responses()
This class is called by responses to make sure that all responses requested are implemented in either...
 
virtual void global_responses(std::vector< double > &resp, const FuelCell::ApplicationCore::FEVector &sol)
This class is used to evaluate all responses that do not require looping over cells. 
 
~OptimizationBlockMatrixApplication()
Destructor. 
Definition: optimization_block_matrix_application.h:92
 
virtual void declare_parameters(ParameterHandler ¶m)
Declare all parameters that are needed for: 
 
unsigned int get_n_dvar() const 
Member function that returns the number of design variables. 
 
virtual void cell_responses(std::vector< double > &resp, const typename DoFApplication< dim >::CellInfo &info, const FuelCell::ApplicationCore::FEVector &sol)
This class is used to evaluate all the functionals that require looping over cells. 
 
const std::vector< std::string > get_all_responses_names() const 
Function used to return all possible response names in OpenFCST. 
Definition: optimization_block_matrix_application.h:293
 
Objects of this kind are used to notify interior applications of changes provoked by an outer loop...
Definition: event.h:51
 
std::vector< std::string > get_name_responses() const 
Member function that returns the name of the responses. 
Definition: optimization_block_matrix_application.h:284
 
unsigned int n_dvar
Number of design variables. 
Definition: optimization_block_matrix_application.h:436
 
std::vector< std::string > get_name_dvar() const 
Member function that returns the name of the design variables. 
Definition: optimization_block_matrix_application.h:275
 
BlockVector< double > FEVector
The vector class used by applications. 
Definition: application_data.h:46
 
The data type used in function calls of Application. 
Definition: fe_vectors.h:59
 
void set_all_response_names()
Initialize the data member all_response_names with all the available Responses in OpenFCST...
 
Application handling matrices and assembling the linear system to solve the sensitivity equations...
Definition: optimization_block_matrix_application.h:49
 
bool optimization
Decision variable for whether the application is being used for optimization. 
Definition: optimization_block_matrix_application.h:465
 
void print_dresponses_dl(std::vector< std::vector< double > > pdf_pdl)
Auxiliary routine to print the values df_dl This routine should be called once df_df is assembled...
Definition: optimization_block_matrix_application.h:366
 
Base class for all linear applications, i.e., all applications requiring Triangulation and DoFHandler...
Definition: dof_application.h:107
 
boost::shared_ptr< ApplicationData > data
Object for auxiliary data. 
Definition: application_base.h:348
 
virtual void cell_dresponses_dl(std::vector< std::vector< double > > &cell_df_dl, const typename DoFApplication< dim >::CellInfo &info, const FuelCell::ApplicationCore::FEVector &src)
This class is used to evaluate the derivative of all the functionals that require looping over cells ...
 
static const FuelCell::ApplicationCore::Event sensitivity_analysis
The Event set by OptimizationMatrixApplication if if we want to compute the sensitivities and a new m...
Definition: optimization_block_matrix_application.h:59
 
const std::string extend_filename(const std::string &, const int precision=3) const 
Member function used in order to extend the name of a file with the design variable name and value...
 
unsigned int n_resp
Number of responses = n_obj + n_con. 
Definition: optimization_block_matrix_application.h:444
 
virtual void bdry_responses(std::vector< double > &resp, const typename DoFApplication< dim >::FaceInfo &info, const FuelCell::ApplicationCore::FEVector &sol)
This class is used to evaluate all the functionals that require looping over boundaries.