16 #ifndef _FUELCELL__ADAPTIVEREFINEMENT_H_
17 #define _FUELCELL__ADAPTIVEREFINEMENT_H_
20 #include <deal.II/base/parameter_handler.h>
21 #include <deal.II/base/convergence_table.h>
22 #include <deal.II/base/smartpointer.h>
23 #include <deal.II/grid/tria.h>
38 namespace ApplicationCore
116 void solve(
const std::string param_file,
117 ParameterHandler& param,
118 bool solution_component_changes_between_data_files =
false );
123 void run_app(
bool solution_component_changes_between_data_files =
false);
128 void run_app ( std::vector<double>& resp,
129 bool solution_component_changes_between_data_files =
false );
134 void run_app ( std::vector<double>& resp,
135 std::vector<std::vector<double> >& dresp_dl,
136 bool solution_component_changes_between_data_files =
false );
143 const std::string dvar,
145 std::vector<double>& resp,
146 std::vector<std::vector<double> >& dresp,
147 const bool gradient =
true );
bool output_initial_mesh
Flag to specify if the initial grid should be saved to a file.
Definition: adaptive_refinement.h:198
unsigned int n_ref
Number of initial refinements for the original mesh.
Definition: adaptive_refinement.h:269
bool L1_L2_error_and_convergence_rate
Use true if you have the exact or analytical solution to compare your numerical results with...
Definition: adaptive_refinement.h:252
std::string filename_initial_mesh
Filename where to output the initial grid.
Definition: adaptive_refinement.h:188
bool nonlinear_solver_for_linear_problem
Use true along with zero initial guess defined in your application to solve a linear problem by means...
Definition: adaptive_refinement.h:262
const FuelCell::ApplicationCore::FEVector & get_coarse_solution() const
This function returns coarse_solution.
Definition: adaptive_refinement.h:167
bool output_intermediate_sol
Boolean flag that is set to true if intermediate solutions, i.e.
Definition: adaptive_refinement.h:212
FuelCell::ApplicationCore::FEVector coarse_solution
Global FE solution at the support points of the initial mesh.
Definition: adaptive_refinement.h:299
void declare_parameters(ParameterHandler ¶m) const
Declare all parameters that are needed for:
FuelCell::ApplicationCore::FEVector solution
Global FE solution at the support points of the most refined mesh.
Definition: adaptive_refinement.h:294
This class is initialized with an application that describes the linearization of the problem that we...
Definition: adaptive_refinement.h:52
FuelCell::ApplicationCore::OptimizationBlockMatrixApplication< dim > * app_linear
Pointer to linear application.
Definition: adaptive_refinement.h:284
FuelCell::ApplicationCore::ApplicationWrapper * app
Poiner to nonlinear application.
Definition: adaptive_refinement.h:289
AdaptiveRefinement(FuelCell::ApplicationCore::ApplicationBase &application)
Constructor.
Definition: adaptive_refinement.h:60
void print_convergence_table()
const FuelCell::ApplicationCore::FEVector & get_solution() const
This function returns solution.
Definition: adaptive_refinement.h:158
void initialize(ParameterHandler ¶m)
Set up how many equations are needed and read in parameters for the parameter handler in order to ini...
std::vector< ConvergenceTable > convergence_tables
If the exact or analytical solution is available, then those convergence tables collect all necessary...
Definition: adaptive_refinement.h:307
Base class for applications.
Definition: application_base.h:113
This class implements either iterative or time-stepping wrapper of applications.
Definition: application_wrapper.h:39
bool output_final_sol
Boolean flag used to specify if the final solution should be output to a file.
Definition: adaptive_refinement.h:223
void solve(const std::string param_file, ParameterHandler ¶m, bool solution_component_changes_between_data_files=false)
Solve the nonlinear problem.
bool gradients
Compute the gradients?
Definition: adaptive_refinement.h:274
std::string filename_final_sol
Filename where the final solution will be output.
Definition: adaptive_refinement.h:231
bool output_intermediate_resp
Bool flag used to specify if all responses/functionals should be evaluated at intermediate refinement...
Definition: adaptive_refinement.h:245
const Triangulation< dim > & get_coarse_triangulation() const
This function returns coarse_triangulation.
Definition: adaptive_refinement.h:176
void run_app(bool solution_component_changes_between_data_files=false)
Run application.
Triangulation< dim > coarse_triangulation
Initial mesh.
Definition: adaptive_refinement.h:279
void test_derivatives(const std::string input_file, const std::string dvar, const double value, std::vector< double > &resp, std::vector< std::vector< double > > &dresp, const bool gradient=true)
Member function used to test the derivatives.
BlockVector< double > FEVector
The vector class used by applications.
Definition: application_data.h:46
Application handling matrices and assembling the linear system to solve the sensitivity equations...
Definition: optimization_block_matrix_application.h:49
void print_parameters() const
Print parameters: