OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | List of all members
AppFrame::OptimizationBlockMatrixApplication< dim > Class Template Referenceabstract

Application handling matrices and assembling the linear system to solve the sensitivity equations. More...

#include <optimization_block_matrix_application.h>

Inheritance diagram for AppFrame::OptimizationBlockMatrixApplication< dim >:
Inheritance graph
[legend]
Collaboration diagram for AppFrame::OptimizationBlockMatrixApplication< dim >:
Collaboration graph
[legend]

Public Member Functions

 OptimizationBlockMatrixApplication (boost::shared_ptr< AppFrame::ApplicationData > data=boost::shared_ptr< AppFrame::ApplicationData >())
 Constructor for an object, owning its own mesh and dof handler.
 
 OptimizationBlockMatrixApplication (AppFrame::DoFApplication< dim > &, bool triangulation_only)
 Constructor for an object, borrowing mesh and dof handler from another object.
 
 ~OptimizationBlockMatrixApplication ()
 Destructor.
 
virtual void declare_parameters (ParameterHandler &param)
 Declare all parameters that are needed for:
 
virtual void set_parameters (const std::vector< std::string > &name_dvar, const std::vector< double > &value_dvar, ParameterHandler &param)
 Function called by optimization loop in order to set the values in the ParameterHandler to the new design parameters.
 
virtual void initialize (ParameterHandler &param)
 Initialize application.
 
virtual void init_solution (AppFrame::FEVector &src)=0
 Initialize nonlinear solution.
 
virtual void responses (std::vector< double > &f, const AppFrame::FEVectors &vectors)
 Post-processing.
 
virtual void check_responses ()
 This class is called by responses to make sure that all responses requested are implemented in either cell_responses, global_responses or face_responses.
 
virtual void cell_responses (std::vector< double > &resp, const typename DoFApplication< dim >::CellInfo &info, const AppFrame::FEVector &sol)
 This class is used to evaluate all the functionals that require looping over cells.
 
virtual void global_responses (std::vector< double > &resp, const AppFrame::FEVector &sol)
 This class is used to evaluate all responses that do not require looping over cells.
 
virtual void print_responses (std::vector< double > &resp)
 This function is used to print the responses to screen (deallog)
 
virtual void dresponses_dl (std::vector< std::vector< double > > &df_dl, const AppFrame::FEVectors &src)
 Post-processing.
 
virtual void cell_dresponses_dl (std::vector< std::vector< double > > &cell_df_dl, const typename DoFApplication< dim >::CellInfo &info, const AppFrame::FEVector &src)
 This class is used to evaluate the derivative of all the functionals that require looping over cells with respect to the design variables.
 
virtual void global_dresponses_dl (std::vector< std::vector< double > > &df_dl, const AppFrame::FEVector &src)
 This class is used to evaluate the sensitivities of all responses that do not require looping over cells with respect of the design variables.
 
virtual void dresponses_du (std::vector< AppFrame::FEVector > &dst, const AppFrame::FEVectors &src)
 Loop over all cells and assemble the vector

\[ \frac{\partial f_m}{\partial u_i} \]

that is used to solve the sensitivity equations by using the local matrix integration functions cell_dfunctional_du(), bdry_dfunctional_dlu() and edge_dfunctinal_du() provided by the derived class.

 
virtual void cell_dresponses_du (std::vector< AppFrame::FEVector > &df_du, const typename DoFApplication< dim >::CellInfo &info, std::vector< std::vector< double > > &src)
 Integration of local bilinear form.
 
virtual void global_dresponses_du (std::vector< AppFrame::FEVector > &df_du, const AppFrame::FEVector &src)
 This class is used to evaluate the sensitivities of all responses that do not require looping over cells with respect of the design variables.
 
virtual void dresidual_dlambda (std::vector< AppFrame::FEVector > &dst, const AppFrame::FEVectors &src)
 Loop over all cells and assemble the vector

\[ \frac{\partial R_n}{\partial \lambda_k} \]

that is used to solve the sensitivity equations by using the local matrix integration functions cell_dresidual_dlambda(), bdry_dresidual_dlambda() and edge_dresidual_dlambda() provided by the derived class.

 
virtual void cell_dresidual_dlambda (std::vector< AppFrame::FEVector > &cell_vector, const typename DoFApplication< dim >::CellInfo &cell, std::vector< std::vector< double > > &src)
 Integration of local bilinear form.
 
void solve_direct (std::vector< std::vector< double > > &df_dl, const AppFrame::FEVectors &sol)
 Solver in order to obtain the analytical sensitivities for the given objective and constraints using the direct method.
 
void solve_adjoint (std::vector< std::vector< double > > &df_dl, const AppFrame::FEVector &sol)
 Solver in order to obtain the analytical sensitivities for the given objective and constraints using the adjoint method.
 
unsigned int get_n_resp () const
 Member function that returns the number of responses.
 
unsigned int get_n_dvar () const
 Member function that returns the number of design variables.
 
std::vector< std::string > get_name_dvar () const
 Member function that returns the name of the design variables.
 
std::vector< std::string > get_name_responses () const
 Member function that returns the name of the responses.
 
void print_default_parameter_file ()
 Print the default parameter handler file.
 
bool get_bool_transfer_solution ()
 Function to see if we are transferring a solution on a refined grid to the initial coarse grid.
 
virtual void set_read_initial_solution (bool &read)
 
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 variables n_dvar, n_resp, name_design_var and name_response.
 
void set_output_variables (std::vector< std::string > &dakota_name_responses)
 
- Public Member Functions inherited from AppFrame::BlockMatrixApplication< dim >
void do_assemble (const FEVectors &, typename DoFHandler< dim >::active_cell_iterator begin, typename DoFHandler< dim >::active_cell_iterator end, Threads::Mutex &mutex)
 Multithreaded assemble function called by assemble().
 
virtual void post_cell_assemble ()
 A call back function for assemble(), called after all cell matrices have been entered, but before the face matrices are computed.
 
 BlockMatrixApplication (boost::shared_ptr< ApplicationData > data=boost::shared_ptr< ApplicationData >())
 Constructor for an object, owning its own mesh and dof handler.
 
 BlockMatrixApplication (DoFApplication< dim > &, bool triangulation_only)
 Constructor for an object, borrowing mesh and dof handler from another object.
 
void _initialize (ParameterHandler &param)
 Initialize data of this class.
 
void add_vector_for_cell_matrix (std::string name, unsigned int block, unsigned int components, bool values, bool derivatives)
 Add a named vector from data to be handed to the cell_matrix() function via CellInfo::values and CellInfo::derivatives.
 
void add_vector_for_flux_matrix (std::string name, unsigned int block, unsigned int components, bool values, bool derivatives)
 Add a named vector from data to be handed to the bdry_matrix() and face_matrix() functions via FaceInfo::values and FaceInfo::derivatives.
 
void remesh_matrices ()
 Initialize sparsity patterns and matrices for the new mesh.
 
virtual void remesh ()
 Refine grid accordingly to the Refinement options under "Grid Generation" in the parameter file.
 
void assemble (const FEVectors &)
 Loop over all cells and assemble the system #matrices by using the local matrix integration functions cell_matrix(), bdry_matrix() and face_matrix() provided by the derived class.
 
void assemble_numerically (const FEVectors &src, const double delta=1e-6)
 Compute the Jacobian of the system of equations,

\[ J(i,j) = \frac{\partial R_i}{\partial u_j} \]

by using forward differences.

 
void residual_constraints (FEVector &dst) const
 Redefinition of residual_constraints() in DoFHandler.
 
virtual void cell_matrix (MatrixVector &cell_matrices, const typename DoFApplication< dim >::CellInfo &cell)
 Integration of local bilinear form.
 
virtual void bdry_matrix (MatrixVector &face_matrices, const typename DoFApplication< dim >::FaceInfo &face)
 Integration of local bilinear form.
 
virtual void face_matrix (MatrixVector &matrices11, MatrixVector &matrices12, MatrixVector &matrices21, MatrixVector &matrices22, const typename DoFApplication< dim >::FaceInfo &face1, const typename DoFApplication< dim >::FaceInfo &face2)
 Integration of local bilinear form.
 
virtual void dirichlet_bc (std::map< unsigned int, double > &boundary_values) const
 
- Public Member Functions inherited from AppFrame::DoFApplication< dim >
 DoFApplication (boost::shared_ptr< ApplicationData > data=boost::shared_ptr< ApplicationData >(), bool multigrid=false)
 Constructor for an object, owning its own mesh and dof handler.
 
 DoFApplication (bool multigrid)
 Constructor for an object owning its own mesh and dof handler and creating new ApplicationData.
 
 DoFApplication (DoFApplication< dim > &, bool triangulation_only, bool multigrid=false)
 Constructor for an object, borrowing mesh and dof handler from another object.
 
 ~DoFApplication ()
 Destructor which deletes owned objects.
 
void _initialize (ParameterHandler &param)
 Initialize from parameter values.
 
void add_vector_for_transfer (FEVector *)
 Add the vector to be transfered from one mesh to the next.
 
void delete_vector_for_transfer ()
 Delete the vector to be transfered from one mesh to the next.
 
virtual void remesh_dofs ()
 Initialize dof handler, count the dofs in each block and renumber the dofs.
 
virtual void init_vector (FEVector &dst) const
 Reinitialize the BlockVector dst such that it contains block_size.size() blocks.
 
virtual double evaluate (const FEVectors &)
 Evaluate a functional during postprocessing such as drag or lift on an aerodynamics problem.
 
virtual double estimate (const FEVectors &src)
 Estimate the error.
 
virtual double residual (FEVector &dst, const FEVectors &src, bool apply_boundaries=true)
 Loop over all cells to compute the residual.
 
unsigned int memory_consumption () const
 Compute the amount of memory needed by this object.
 
virtual double cell_estimate (const CellInfo &)
 Integration of estimate inside a cell.
 
virtual double bdry_estimate (const FaceInfo &)
 Integration of estimate on a boundary face.
 
virtual double face_estimate (const FaceInfo &, const FaceInfo &)
 Integration of estimate on an interior face.
 
virtual void cell_residual (FEVector &cell_vector, const CellInfo &cell)
 Integration of local linear form.
 
virtual void bdry_residual (FEVector &face_vector, const FaceInfo &face)
 Integration of local linear form.
 
virtual void face_residual (FEVector &face_vector1, FEVector &face_vector2, const FaceInfo &face1, const FaceInfo &face2)
 Integration of local linear form.
 
template<class BOX , class WORKER >
void integrate_1form (BOX &box, WORKER &local_forms, FEVector &dst, const FEVectors &src)
 The function integrating 1-forms as for instance used by residual().
 
void integrate_functional (LocalEstimate< dim > &local_forms, FEVectors &dst, const FEVectors &src)
 The function integrating functionals using local integrators.
 
void store_triangulation (Triangulation< dim > &new_tr)
 Function to copy a triangulation object for use after refinement.
 
void transfer_solution_to_coarse_mesh (Triangulation< dim > &tr_coarse, AppFrame::FEVector &coarse_solution, AppFrame::FEVector &refined_solution)
 Function to perform the transfer of a solution on a refined grid to the initial coarse grid.
 
void read_init_solution (AppFrame::FEVector &dst, bool &good_solution) const
 Function to transfer a solution on a refined grid to the initial coarse grid.
 
template<class INFO , typename TYPE >
void fill_local_data (const INFO &info, const std::vector< VectorSelector > &data_vector, std::vector< std::vector< std::vector< TYPE > > > &data) const
 Fill local data vector with values of the stored finite element function vectors.
 
virtual void grid_out (const std::string &basename)
 Output the grid used to solve the problem.
 
virtual void data_out (const std::string &basename, const FEVectors &src)
 
virtual void data_out (const std::string &basename, const FEVectors &src, const std::vector< std::string > &solution_names)
 
virtual void data_out (const std::string &basename, const FEVector &solution, const std::vector< std::string > &solution_names, const FEVector &postprocessing=FEVector(), const std::vector< std::string > &postprocessing_names=std::vector< std::string >())
 This function outputs the results of a computation.
 
- Public Member Functions inherited from AppFrame::ApplicationBase
 ApplicationBase (boost::shared_ptr< ApplicationData > ex_data=boost::shared_ptr< ApplicationData >())
 Constructor for an application.
 
 ApplicationBase (const ApplicationBase &other)
 Copy constructor.
 
virtual ~ApplicationBase ()
 Virtual destructor.
 
void print_parameters_to_file (ParameterHandler &param, const std::string &file_name, const ParameterHandler::OutputStyle &style)
 Print default parameters for the application to a file.
 
virtual void start_vector (FEVector &dst, std::string caller) const
 Initialize vector to problem size.
 
virtual void solve (FEVector &start, const FEVectors &rhs)
 Solve the system assembled with right hand side rhs and return the result in start.
 
virtual void Tsolve (FEVector &start, const FEVectors &rhs)
 Solve the dual system assembled with right hand side rhs and return the result in start.
 
virtual void grid_out (const std::string &filename) const
 Write the mesh in the format specified by the ParameterHandler.
 
virtual void data_out (const std::string &filename, const FEVectors &src, const std::vector< std::string >)
 
boost::shared_ptr
< ApplicationData
get_data ()
 Get access to the protected variable data.
 
const boost::shared_ptr
< ApplicationData
get_data () const
 Get read-only access to the protected variable data.
 
virtual std::string id () const
 Return a unique identification string for this application.
 
virtual void notify (const Event &reason)
 Add a reason for assembling.
 
virtual void clear ()
 Reset the application class such that a call to initialize() is possible again and will produce the same result as if the object was fresh.
 
void clear_events ()
 Clear all notifications.
 

Static Public Attributes

static const AppFrame::Event sensitivity_analysis
 The Event set by OptimizationMatrixApplication if if we want to compute the sensitivities and a new matrix should be assembled.
 

Protected Types

typedef void(* CELL_Dvalues )(std::vector< AppFrame::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 of either residual or responses with respect to design variables.
 

Protected Member Functions

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.
 
void print_dresponses_du (std::vector< AppFrame::FEVector > df_du)
 Auxiliary routine to print the values of df_du This routine should be called once df_du is assembled.
 
void dfunction (std::vector< AppFrame::FEVector > &dst, const AppFrame::FEVectors &src, bool dfunctional_du, bool dresidual_dlambda)
 This is an auxiliary function called by dresidual_dlambda and dfunctional_du.
 

Protected Attributes

unsigned int n_dvar
 Number of design variables.
 
unsigned int n_obj
 Number of objective functions.
 
unsigned int n_resp
 Number of responses = n_obj + n_con.
 
std::vector< std::string > name_design_var
 Member that stores the name of the design variables.
 
std::vector< std::string > name_responses
 Member that stores the name of the responses, i.e.
 
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 optimization is not being used and ignored otherwise.
 
bool output_coarse_solution
 Decision variable for whether the solution is to be output for transfer.
 
bool read_in_initial_solution
 Decision variable for whether the initial solution should be read from file.
 
bool optimization
 Decision variable for whether the application is being used for optimization.
 
- Protected Attributes inherited from AppFrame::BlockMatrixApplication< dim >
Table< 2, DoFTools::Coupling > cell_couplings
 Couplings through cell bilinear forms.
 
Table< 2, DoFTools::Coupling > flux_couplings
 Couplings through flux bilinear forms.
 
boost::shared_ptr< Quadrature
< dim > > 
quadrature_assemble_cell
 Quadrature rule for matrix assembling on cells.
 
boost::shared_ptr< Quadrature
< dim-1 > > 
quadrature_assemble_face
 Quadrature rule for matrix assembling on faces.
 
BlockSparseMatrix< double > matrix
 Storage for the actual matrix.
 
FilteredMatrix< BlockVector
< double > > 
filtered_matrix
 Wrapper for the actual matrix that contains all the boundary conditions.
 
- Protected Attributes inherited from AppFrame::DoFApplication< dim >
std::map< unsigned int, double > boundary_constraints
 List of all nodes constrained by a strong boundary condition, together with a value to be assigned.
 
ConstraintMatrix hanging_node_constraints
 Constraint Matrix object.
 
boost::shared_ptr< Mapping< dim > > mapping
 The mapping used for the transformation of mesh cells.
 
unsigned int mapping_degree
 Degree used for polynomial mapping of arbitrary order.
 
boost::shared_ptr
< Triangulation< dim > > 
tr
 Pointer to the Triangulation object.
 
boost::shared_ptr
< FiniteElement< dim > > 
element
 The finite element used in dof.
 
boost::shared_ptr< DoFHandler
< dim > > 
dof
 Pointer to the DoFHandler object.
 
boost::shared_ptr
< MGDoFHandler< dim > > 
mg_dof
 Pointer to the MGDoFHandler object.
 
MeshWorker::BlockInfo block_info
 
Vector< float > cell_errors
 The result of error estimation by cell.
 
Vector< float > face_errors
 The result of error estimation by face.
 
GridOut g_out
 The object for writing grids.
 
DataOut< dim, DoFHandler< dim > > d_out
 The object for writing data.
 
std::vector
< DataComponentInterpretation::DataComponentInterpretation > 
data_interpretation
 Vector for adjusting output to vector data.
 
std::string refinement
 Refinement parameter from parameter file.
 
unsigned int initial_refinement
 Initial refinement from parameter file.
 
double refinement_threshold
 Refinement threshold for adaptive method from parameter file.
 
double coarsening_threshold
 Coarsening threshold for adaptive method from parameter file.
 
bool sort_cuthill
 Flag for sorting with Cuthill McKee algorithm.
 
Point< dimsort_direction
 Direction for downstream sorting.
 
std::vector< FEVector * > transfer_vectors
 List of vector names to be transfered from one grid to the next.
 
Quadrature< dimquadrature_residual_cell
 Quadrature rule for residual computation on cells.
 
Quadrature< dim-1 > quadrature_residual_bdry
 Quadrature rule for residual computation on boundary faces.
 
Quadrature< dim-1 > quadrature_residual_face
 Quadrature rule for residual computation on faces.
 
bool boundary_fluxes
 Extend the integration loops in assemble() and residual() also to boundary faces.
 
bool interior_fluxes
 Extend the integration loops in assemble() and residual() also to interior faces.
 
- Protected Attributes inherited from AppFrame::ApplicationBase
boost::shared_ptr
< ApplicationData
data
 Object for auxiliary data.
 
Event notifications
 Accumulate reasons for reassembling here.
 

Additional Inherited Members

- Public Types inherited from AppFrame::DoFApplication< dim >
typedef
MeshWorker::InfoObjects::IntegrationInfo
< dim, FEValuesBase< dim > > 
CellInfo
 This is a redefinition of CellInfo in order to call MeshWorker objects instead of the previous implementation in AppFrame.
 
typedef
MeshWorker::InfoObjects::IntegrationInfo
< dim, FEFaceValuesBase< dim > > 
FaceInfo
 This is a redefinition of CellInfo in order to call MeshWorker objects instead of the previous implementation in AppFrame.
 
- Public Attributes inherited from AppFrame::DoFApplication< dim >
unsigned int verbosity
 Controls verbosity of certain functions.
 
std::vector
< DataComponentInterpretation::DataComponentInterpretation > 
solution_interpretations
 solution_interpretations identifies whether some solution_names are scalars or parts of a vector.
 
std::vector
< DataComponentInterpretation::DataComponentInterpretation > 
postprocessing_interpretations
 postprocessing_interpretations identifies whether some postprocessing_names are scalars or parts of a vector.
 
bool output_materials_and_levels
 output_materials_and_levels if true then visualized, otherwise suppressed.
 
Vector< double > output_materials
 Vector that will be used by data_out to store the material ids.
 
Vector< double > output_levels
 Vector that will be used by data_out to store the refinement levels at each cell.
 

Detailed Description

template<int dim>
class AppFrame::OptimizationBlockMatrixApplication< dim >

Application handling matrices and assembling the linear system to solve the sensitivity equations.

This class inherits the information on Triangulation and DoFHandler from its base class BlockMatrixApplication. It adds information on the linear system for sensitivity analysis

Author
Marc Secanell, 2006

Member Typedef Documentation

template<int dim>
typedef void(* AppFrame::OptimizationBlockMatrixApplication< dim >::CELL_Dvalues)(std::vector< AppFrame::FEVector > &, const typename DoFApplication< dim >::CellInfo &, std::vector< std::vector< double > > &)
protected

Definition of a pointer to a function that is called in a loop over cells to compute the derivatives of either residual or responses with respect to design variables.

Constructor & Destructor Documentation

Constructor for an object, owning its own mesh and dof handler.

If data is null, a new handler for application data is generated.

see constructor of DoFApplication.

Constructor for an object, borrowing mesh and dof handler from another object.

If triangulation_only is true, only the triangulation is borrowed.

see constructor of DoFApplication.

Destructor.

Member Function Documentation

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::cell_dresidual_dlambda ( std::vector< AppFrame::FEVector > &  cell_vector,
const typename DoFApplication< dim >::CellInfo cell,
std::vector< std::vector< double > > &  src 
)
virtual

Integration of local bilinear form.

Reimplemented in FuelCell::Application::AppCathode< dim >.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::cell_dresponses_dl ( std::vector< std::vector< double > > &  cell_df_dl,
const typename DoFApplication< dim >::CellInfo info,
const AppFrame::FEVector src 
)
virtual

This class is used to evaluate the derivative of all the functionals that require looping over cells with respect to the design variables.

This class is called by responses to evaluate the response at each cell.

Reimplemented in FuelCell::Application::AppCathode< dim >, and FuelCell::Application::AppPemfc< dim >.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::cell_dresponses_du ( std::vector< AppFrame::FEVector > &  df_du,
const typename DoFApplication< dim >::CellInfo info,
std::vector< std::vector< double > > &  src 
)
virtual

Integration of local bilinear form.

Reimplemented in FuelCell::Application::AppCathode< dim >, and FuelCell::Application::AppPemfc< dim >.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::cell_responses ( std::vector< double > &  resp,
const typename DoFApplication< dim >::CellInfo info,
const AppFrame::FEVector sol 
)
virtual

This class is used to evaluate all the functionals that require looping over cells.

This class is called by responses to evaluate the response at each cell.

Reimplemented in FuelCell::Application::AppCathode< dim >, FuelCell::Application::AppPemfc< dim >, FuelCell::Application::AppLaplace< dim >, and FuelCell::Application::AppReadMesh< dim >.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::check_responses ( )
virtual

This class is called by responses to make sure that all responses requested are implemented in either cell_responses, global_responses or face_responses.

Reimplemented in FuelCell::Application::AppCathode< dim >, and FuelCell::Application::AppPemfc< dim >.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::declare_parameters ( ParameterHandler &  param)
virtual

Declare all parameters that are needed for:

  • the computation of the equation coefficients
  • the control of the linear system solution
  • ...

Reimplemented from AppFrame::BlockMatrixApplication< dim >.

Reimplemented in FuelCell::Application::AppCathode< dim >, FuelCell::Application::AppPemfc< dim >, FuelCell::Application::AppReadMesh< dim >, and FuelCell::Application::AppLaplace< dim >.

Referenced by AppFrame::OptimizationBlockMatrixApplication< dim >::print_default_parameter_file().

Here is the caller graph for this function:

template<int dim>
void AppFrame::OptimizationBlockMatrixApplication< dim >::dfunction ( std::vector< AppFrame::FEVector > &  dst,
const AppFrame::FEVectors src,
bool  dfunctional_du,
bool  dresidual_dlambda 
)
protected

This is an auxiliary function called by dresidual_dlambda and dfunctional_du.

Since both functions have to do the same thing, loop over the cells, but they call a different member function for cell assembly, I use a flag to change from one to the other and then I call this generic function. With this dresidual_dlambda and dfunctional_du have the same interface and the changes happen to a protected class.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::dresidual_dlambda ( std::vector< AppFrame::FEVector > &  dst,
const AppFrame::FEVectors src 
)
virtual

Loop over all cells and assemble the vector

\[ \frac{\partial R_n}{\partial \lambda_k} \]

that is used to solve the sensitivity equations by using the local matrix integration functions cell_dresidual_dlambda(), bdry_dresidual_dlambda() and edge_dresidual_dlambda() provided by the derived class.

The vector index refers to k and the FEVector index refers to n.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::dresponses_dl ( std::vector< std::vector< double > > &  df_dl,
const AppFrame::FEVectors src 
)
virtual

Post-processing.

Evaluate the partial derivative of all functionals necessary such as the objective function and the constraints of the optimization problem In this member class we obtain the term:

\[ \frac{\partial f_m}{\partial \lambda_k} \]

Note
size of df_dl is for f num_objectives + num_constraints and for lambda = num_design_var
template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::dresponses_du ( std::vector< AppFrame::FEVector > &  dst,
const AppFrame::FEVectors src 
)
virtual

Loop over all cells and assemble the vector

\[ \frac{\partial f_m}{\partial u_i} \]

that is used to solve the sensitivity equations by using the local matrix integration functions cell_dfunctional_du(), bdry_dfunctional_dlu() and edge_dfunctinal_du() provided by the derived class.

The vector index refers to m and the BlockVector index refers to i.

template<int dim>
bool AppFrame::OptimizationBlockMatrixApplication< dim >::get_bool_transfer_solution ( )
inline

Function to see if we are transferring a solution on a refined grid to the initial coarse grid.

References AppFrame::OptimizationBlockMatrixApplication< dim >::output_coarse_solution.

template<int dim>
unsigned int AppFrame::OptimizationBlockMatrixApplication< dim >::get_n_dvar ( ) const

Member function that returns the number of design variables.

template<int dim>
unsigned int AppFrame::OptimizationBlockMatrixApplication< dim >::get_n_resp ( ) const

Member function that returns the number of responses.

template<int dim>
std::vector<std::string> AppFrame::OptimizationBlockMatrixApplication< dim >::get_name_dvar ( ) const
inline

Member function that returns the name of the design variables.

References AppFrame::OptimizationBlockMatrixApplication< dim >::name_design_var.

template<int dim>
std::vector<std::string> AppFrame::OptimizationBlockMatrixApplication< dim >::get_name_responses ( ) const
inline

Member function that returns the name of the responses.

References AppFrame::OptimizationBlockMatrixApplication< dim >::name_responses.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::global_dresponses_dl ( std::vector< std::vector< double > > &  df_dl,
const AppFrame::FEVector src 
)
virtual

This class is used to evaluate the sensitivities of all responses that do not require looping over cells with respect of the design variables.

An example of one of this types of constraints is the solid volume fraction.

Reimplemented in FuelCell::Application::AppCathode< dim >, and FuelCell::Application::AppPemfc< dim >.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::global_dresponses_du ( std::vector< AppFrame::FEVector > &  df_du,
const AppFrame::FEVector src 
)
virtual

This class is used to evaluate the sensitivities of all responses that do not require looping over cells with respect of the design variables.

An example of one of this types of constraints is the solid volume fraction.

Reimplemented in FuelCell::Application::AppCathode< dim >, and FuelCell::Application::AppPemfc< dim >.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::global_responses ( std::vector< double > &  resp,
const AppFrame::FEVector sol 
)
virtual

This class is used to evaluate all responses that do not require looping over cells.

An example of one of this types of constraints is the solid volume fraction.

Reimplemented in FuelCell::Application::AppCathode< dim >, FuelCell::Application::AppPemfc< dim >, FuelCell::Application::AppLaplace< dim >, and FuelCell::Application::AppReadMesh< dim >.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::init_solution ( AppFrame::FEVector src)
pure virtual
template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::initialize ( ParameterHandler &  param)
virtual
template<int dim>
void AppFrame::OptimizationBlockMatrixApplication< dim >::print_default_parameter_file ( )
inline

Print the default parameter handler file.

References AppFrame::OptimizationBlockMatrixApplication< dim >::declare_parameters().

Here is the call graph for this function:

template<int dim>
void AppFrame::OptimizationBlockMatrixApplication< dim >::print_dresponses_dl ( std::vector< std::vector< double > >  pdf_pdl)
inlineprotected

Auxiliary routine to print the values df_dl This routine should be called once df_df is assembled.

This is done on solve_direct

References AppFrame::OptimizationBlockMatrixApplication< dim >::n_dvar, and AppFrame::OptimizationBlockMatrixApplication< dim >::n_resp.

template<int dim>
void AppFrame::OptimizationBlockMatrixApplication< dim >::print_dresponses_du ( std::vector< AppFrame::FEVector df_du)
inlineprotected

Auxiliary routine to print the values of df_du This routine should be called once df_du is assembled.

This is done on solve_direct

References AppFrame::OptimizationBlockMatrixApplication< dim >::n_resp.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::print_responses ( std::vector< double > &  resp)
virtual

This function is used to print the responses to screen (deallog)

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::responses ( std::vector< double > &  f,
const AppFrame::FEVectors vectors 
)
virtual

Post-processing.

Evaluate all responses, usually functionals, necessary such as the objective function and the constraints of the optimization problem

Note
size of f is num_objectives + num_constraints
template<int dim>
void AppFrame::OptimizationBlockMatrixApplication< dim >::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 variables n_dvar, n_resp, name_design_var and name_response.

template<int dim>
void AppFrame::OptimizationBlockMatrixApplication< dim >::set_output_variables ( std::vector< std::string > &  dakota_name_responses)
inline
template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::set_parameters ( const std::vector< std::string > &  name_dvar,
const std::vector< double > &  value_dvar,
ParameterHandler &  param 
)
virtual

Function called by optimization loop in order to set the values in the ParameterHandler to the new design parameters.

Since ParameterHandler depends on the problem we are solving, set_parameters() is set at the most inner loop of the application.

Reimplemented in FuelCell::Application::AppCathode< dim >, FuelCell::Application::AppPemfc< dim >, and FuelCell::Application::AppLaplace< dim >.

template<int dim>
virtual void AppFrame::OptimizationBlockMatrixApplication< dim >::set_read_initial_solution ( bool &  read)
inlinevirtual
template<int dim>
void AppFrame::OptimizationBlockMatrixApplication< dim >::solve_adjoint ( std::vector< std::vector< double > > &  df_dl,
const AppFrame::FEVector sol 
)

Solver in order to obtain the analytical sensitivities for the given objective and constraints using the adjoint method.

This member function calls the different assemble member functions and solves the system of equations

Note
The adjoint method is the most computationally efficient method if the number of number of design variables is larger than the number of objectives and constraints
template<int dim>
void AppFrame::OptimizationBlockMatrixApplication< dim >::solve_direct ( std::vector< std::vector< double > > &  df_dl,
const AppFrame::FEVectors sol 
)

Solver in order to obtain the analytical sensitivities for the given objective and constraints using the direct method.

This member function calls the different assemble member functions and solves the system of equations

Note
The direct method is the most computationally efficient method if the number of objectives and constraints is larger than the number of design variables

Member Data Documentation

template<int dim>
unsigned int AppFrame::OptimizationBlockMatrixApplication< dim >::n_dvar
protected
template<int dim>
unsigned int AppFrame::OptimizationBlockMatrixApplication< dim >::n_obj
protected

Number of objective functions.

template<int dim>
unsigned int AppFrame::OptimizationBlockMatrixApplication< dim >::n_resp
protected
template<int dim>
std::vector<std::string> AppFrame::OptimizationBlockMatrixApplication< dim >::name_design_var
protected

Member that stores the name of the design variables.

Referenced by AppFrame::OptimizationBlockMatrixApplication< dim >::get_name_dvar().

template<int dim>
std::vector<std::string> AppFrame::OptimizationBlockMatrixApplication< dim >::name_output_var
protected

Member that stores the name of the output variables, These names will be written to name_responses if optimization is not being used and ignored otherwise.

template<int dim>
std::vector<std::string> AppFrame::OptimizationBlockMatrixApplication< dim >::name_responses
protected

Member that stores the name of the responses, i.e.

objectives and constraints.

Note
first we put the name of all objectives and then all constraints

Referenced by AppFrame::OptimizationBlockMatrixApplication< dim >::get_name_responses(), and AppFrame::OptimizationBlockMatrixApplication< dim >::set_output_variables().

template<int dim>
bool AppFrame::OptimizationBlockMatrixApplication< dim >::optimization
protected

Decision variable for whether the application is being used for optimization.

template<int dim>
bool AppFrame::OptimizationBlockMatrixApplication< dim >::output_coarse_solution
protected

Decision variable for whether the solution is to be output for transfer.

Referenced by AppFrame::OptimizationBlockMatrixApplication< dim >::get_bool_transfer_solution().

template<int dim>
bool AppFrame::OptimizationBlockMatrixApplication< dim >::read_in_initial_solution
protected

Decision variable for whether the initial solution should be read from file.

Referenced by AppFrame::OptimizationBlockMatrixApplication< dim >::set_read_initial_solution().

template<int dim>
const AppFrame::Event AppFrame::OptimizationBlockMatrixApplication< dim >::sensitivity_analysis
static

The Event set by OptimizationMatrixApplication if if we want to compute the sensitivities and a new matrix should be assembled.


The documentation for this class was generated from the following file: