OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
dof_application.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 //
3 // FCST: Fuel Cell Simulation Toolbox
4 //
5 // Copyright (C) 2006-2009 by Guido Kanschat
6 // Copyright (C) 2006-2014 by Energy Systems Design Laboratory, University of Alberta
7 //
8 // This software is distributed under the MIT License
9 // For more information, see the README file in /doc/LICENSE
10 //
11 // - Class: dof_application.h
12 // - Description: This class implements:
13 // - triangulation
14 // - finite element spaces
15 // - dof handler
16 // - distribution of dofs over triangulation
17 // - residual of nonlinear system of equations
18 // - estimations of cell-wise errors for adaptive refinement
19 // - and more
20 // - Developers: Guido Kanschat, Texas A&M University
21 // Valentin N. Zingan, University of Alberta
22 // Marc Secanell, University of Alberta
23 // Peter Dobson, University of Alberta
24 //
25 // ----------------------------------------------------------------------------
26 
27 #ifndef _FUEL_CELL_APPLICATION_CORE_DOF_APPLICATION_H_
28 #define _FUEL_CELL_APPLICATION_CORE_DOF_APPLICATION_H_
29 
30 #include <base/parameter_handler.h>
31 #include <base/data_out_base.h>
32 #include <base/quadrature_lib.h>
33 
34 #include <lac/block_indices.h>
35 #include <lac/constraint_matrix.h>
36 #include <lac/block_vector.h>
37 
38 #include <grid/tria.h>
39 #include <grid/tria_iterator.h>
40 #include <grid/tria_accessor.h>
41 #include <grid/grid_refinement.h>
42 #include <grid/grid_tools.h>
43 #include <grid/grid_out.h>
44 
45 #include <multigrid/mg_tools.h>
46 
47 #include <fe/fe.h>
48 #include <fe/mapping_q1.h>
49 #include <fe/mapping_q.h>
50 #include <fe/fe_tools.h>
51 #include <fe/fe_values.h>
52 
53 #include <dofs/dof_handler.h>
54 #include <dofs/dof_renumbering.h>
55 #include <dofs/dof_tools.h>
56 
57 #include <numerics/data_out.h>
58 #include <numerics/solution_transfer.h>
59 #include <numerics/error_estimator.h>
60 
61 #include "application_wrapper.h"
62 #include "mesh_loop_info_objects.h"
63 #include "fcst_utilities.h"
65 #include <grid/geometry.h>
66 #include <grid/geometries.h>
67 
68 #include <deal.II/lac/petsc_parallel_vector.h>
69 #include <deal.II/lac/petsc_vector.h>
70 
71 #include <iostream>
72 #include <fstream>
73 #include <string>
74 #include <sstream>
75 #include <typeinfo>
76 
77 using namespace dealii;
78 using namespace FuelCell::ApplicationCore;
79 
80 namespace dealii
81 {
82  template<int, int> class Triangulation;
83  template<int, int> class Mapping;
84  template<int, int> class FiniteElement;
85  template<int, int> class DoFHandler;
86 }
87 
88 namespace FuelCell
89 {
90  namespace ApplicationCore
91  {
92 
121  template<int dim>
123  {
124  public:
125 
130 
135 
138 
145  DoFApplication(boost::shared_ptr<ApplicationData> data = boost::shared_ptr<ApplicationData>(),
146  bool multigrid = false);
147 
152  DoFApplication(bool multigrid);
153 
163  bool triangulation_only,
164  bool multigrid = false);
165 
170  ~DoFApplication();
171 
190  virtual void declare_parameters(ParameterHandler& param);
191 
196  virtual void initialize(ParameterHandler& param);
197 
206  virtual void remesh_dofs();
207 
212  virtual void remesh();
213 
226  virtual void init_vector(FEVector& dst) const;
227 
305  virtual void initialize_solution (FEVector& initial_guess,
306  std::shared_ptr<Function<dim> > initial_function = std::shared_ptr<Function<dim> >());
307 
309 
312 
322  virtual double estimate(const FEVectors& src);
323 
327  virtual double evaluate(const FEVectors& src);
328 
336  virtual double residual(FEVector& dst,
337  const FEVectors& src,
338  bool apply_boundaries = true);
339 
353  virtual void residual_constraints(FEVector& dst) const;
354 
355  #ifdef OPENFCST_WITH_PETSC
356  virtual void residual_constraints(PETScWrappers::MPI::Vector& dst, const std::vector<unsigned int>&) const;
357  #endif
358 
362  void store_triangulation(Triangulation<dim>& new_tr);
363 
367  void add_vector_for_transfer(FEVector* src);
368 
372  void delete_vector_for_transfer();
373 
377  void transfer_solution_to_coarse_mesh(Triangulation<dim>& tr_coarse,
378  FEVector& coarse_solution,
379  FEVector& refined_solution);
380 
384  unsigned int memory_consumption() const;
386 
389 
399  std::string filename_initial_sol;
400 
424 
429 
442 
444 
447 
450  virtual void grid_out(const std::string& basename);
451 
455  virtual void data_out(const std::string& basename,
456  const FEVectors& src);
457 
465  void print(const std::string& basename,
466  const FEVector& src,
467  const std::vector<unsigned int>& src_indices = std::vector<unsigned int>()) const;
469 
472 
476  std::vector< DataComponentInterpretation::DataComponentInterpretation > solution_interpretations;
477 
482  std::vector< DataComponentInterpretation::DataComponentInterpretation > postprocessing_interpretations;
483 
487  std::vector<DataComponentInterpretation::DataComponentInterpretation> data_interpretation;
488 
494 
499  Vector<double> output_materials;
500 
505  Vector<double> output_levels;
506 
517 
525 
533 
541  std::vector<unsigned int> solution_printing_indices;
542 
550  std::vector<unsigned int> postprocessing_printing_indices;
551 
559 
561 
565  boost::shared_ptr< Boundary<dim> > curved_boundary;
566 
570  types::boundary_id curved_bdry_id;
571 
572  protected:
575 
580  void _initialize(ParameterHandler& param);
581 
591  virtual void initialize_triangulation(ParameterHandler& param);
592 
596  //virtual void initialize_grid(ParameterHandler& param){};
597 
614  void read_init_solution(FEVector& dst,
615  bool& good_solution) const;
617 
620 
623  virtual void cell_residual(FEVector& cell_vector,
624  const CellInfo& cell);
625 
629  virtual void bdry_residual(FEVector& face_vector,
630  const FaceInfo& face);
631 
635  virtual void face_residual(FEVector& face_vector1,
636  FEVector& face_vector2,
637  const FaceInfo& face1,
638  const FaceInfo& face2);
639 
643  virtual double cell_estimate(const CellInfo& src);
644 
648  virtual double bdry_estimate(const FaceInfo& src);
649 
653  virtual double face_estimate(const FaceInfo& src1,
654  const FaceInfo& src2);
656 
659 
667  std::vector< component_materialID_value_map > component_materialID_value_maps;
668 
677  std::vector< component_boundaryID_value_map > component_boundaryID_value_maps;
678 
682 
685  boost::shared_ptr< FuelCellShop::Geometry::GridBase<dim> > mesh_generator;
686 
691  boost::shared_ptr<Triangulation<dim> > tr;
692 
694 
697  GridOut g_out;
698 
702  DataOut<dim, DoFHandler<dim> > d_out;
703 
753  virtual void data_out(const std::string& basename,
754  const FEVector& solution,
755  const std::vector<std::string>& solution_names,
756  const std::vector< DataPostprocessor<dim>* >& PostProcessing);
757 
772  virtual void data_out(const std::string& basename,
773 
774  const FEVector& solution,
775  const std::vector<std::string>& solution_names,
776 
777  const FEVector& postprocessing = FEVector(),
778  const std::vector<std::string>& postprocessing_names = std::vector<std::string>());
779 
784  void constrain_boundary(FEVector& v,
785  bool homogeneous) const;
786 
792 
800  std::map<unsigned int, double> boundary_constraints;
801 
818  ConstraintMatrix hanging_node_constraints;
819 
830  boost::shared_ptr<Mapping<dim> > mapping;
831 
853  unsigned int mapping_degree;
854 
858  boost::shared_ptr<FiniteElement<dim> > element;
859 
863  boost::shared_ptr<DoFHandler<dim> > dof;
864 
868  boost::shared_ptr<MGDoFHandler<dim> > mg_dof;
869 
871 
875  Vector<float> cell_errors;
876 
880  Vector<float> face_errors;
881 
885  std::string refinement;
886 
890  unsigned int initial_refinement;
891 
896 
901 
908 
915  Point<dim> sort_direction;
916 
920  std::vector<FEVector*> transfer_vectors;
921 
925  Quadrature<dim> quadrature_residual_cell;
926 
931 
937 
967  Table<2, DoFTools::Coupling> cell_couplings;
968 
973  Table<2, DoFTools::Coupling> flux_couplings;
974 
980 
986 
991  unsigned int verbosity;
992 
993  private:
994 
1000  virtual void sort_dofs(DoFHandler<dim>* dof_handler,
1001  MGDoFHandler<dim>* mg_dof_handler = NULL ) const;
1002 
1013  virtual void distribute_face_to_cell_errors();
1014 
1028  virtual double global_from_local_errors() const;
1029  };
1030 
1031  } // ApplicationCore
1032 
1033 } // FuelCell
1034 
1035 #endif
GridOut g_out
The object for writing grids.
Definition: dof_application.h:697
Definition: dof_application.h:82
FuelCell::ApplicationCore::IntegrationInfo< dim, FEFaceValuesBase< dim > > FaceInfo
Shortcut.
Definition: dof_application.h:134
boost::shared_ptr< FuelCellShop::Geometry::GridBase< dim > > mesh_generator
Grid.
Definition: dof_application.h:685
const unsigned int dim
Definition: fcst_constants.h:24
unsigned int mapping_degree
Degree used for polynomial mapping of arbitrary order.
Definition: dof_application.h:853
Vector< double > output_materials
Vector that will be used by data_out to store the material ids.
Definition: dof_application.h:499
DataOut< dim, DoFHandler< dim > > d_out
The object for writing data.
Definition: dof_application.h:702
ConstraintMatrix hanging_node_constraints
Constraint Matrix object.
Definition: dof_application.h:818
bool output_initial_sol
Flag to output the initial solution used to start the solving process.
Definition: dof_application.h:413
Definition: dof_application.h:83
BlockInfo block_info
Definition: dof_application.h:870
std::vector< component_boundaryID_value_map > component_boundaryID_value_maps
Each entry of this std::vector reflects the following structure (see FuelCell::InitialAndBoundaryData...
Definition: dof_application.h:677
A small structure collecting the different BlockIndices of FEVector vectors (for instance, solution) involved in the computations.
Definition: mesh_loop_info_objects.h:182
bool output_coarse_solution
Bool flag used to specify if the final solution should be stored in the coarse mesh in order to be us...
Definition: dof_application.h:441
Quadrature< dim-1 > quadrature_residual_face
Quadrature rule for residual computation on faces.
Definition: dof_application.h:936
Definition: mesh_loop_info_objects.h:41
std::vector< FEVector * > transfer_vectors
List of vector names to be transfered from one grid to the next.
Definition: dof_application.h:920
Vector< double > output_levels
Vector that will be used by data_out to store the refinement levels at each cell. ...
Definition: dof_application.h:505
bool output_actual_degree
true, if you want to output solution and postprocessing data using actual finite element fields Q_n w...
Definition: dof_application.h:516
bool print_solution
true, if you want to print FEVector solution to a text file.
Definition: dof_application.h:524
std::vector< unsigned int > postprocessing_printing_indices
The indices of the FEVector postprocessing to be printed to a text file.
Definition: dof_application.h:550
bool sort_cuthill
Flag for sorting with Cuthill McKee algorithm.
Definition: dof_application.h:907
boost::shared_ptr< MGDoFHandler< dim > > mg_dof
Pointer to the MGDoFHandler object.
Definition: dof_application.h:868
Definition: dof_application.h:84
std::string refinement
Refinement parameter from parameter file.
Definition: dof_application.h:885
std::vector< DataComponentInterpretation::DataComponentInterpretation > data_interpretation
Definition: dof_application.h:487
FuelCell::SystemManagement system_management
This object knows everything about FCST equations, variables, couplings, etc.
Definition: dof_application.h:791
Vector< float > cell_errors
The result of error estimation by cell.
Definition: dof_application.h:875
bool interior_fluxes
Extend the integration loops in assemble() and residual() also to interior faces. ...
Definition: dof_application.h:985
std::vector< component_materialID_value_map > component_materialID_value_maps
Each entry of this std::vector reflects the following structure (see FuelCell::InitialAndBoundaryData...
Definition: dof_application.h:667
std::string filename_initial_sol
Filename where to output the initial grid.
Definition: dof_application.h:399
boost::shared_ptr< Boundary< dim > > curved_boundary
Curved boundary.
Definition: dof_application.h:565
types::boundary_id curved_bdry_id
Curved boundary ID.
Definition: dof_application.h:570
bool use_predefined_solution
Use user pre-defined initial solution.
Definition: dof_application.h:428
FuelCell::ApplicationCore::IntegrationInfo< dim, FEValuesBase< dim > > CellInfo
Shortcut.
Definition: dof_application.h:129
This class is created for the objects handed to the mesh loops.
Definition: mesh_loop_info_objects.h:625
Vector< float > face_errors
The result of error estimation by face.
Definition: dof_application.h:880
double coarsening_threshold
Coarsening threshold for adaptive method from parameter file.
Definition: dof_application.h:900
boost::shared_ptr< Triangulation< dim > > tr
Pointer to the Triangulation object.
Definition: dof_application.h:691
Quadrature< dim > quadrature_residual_cell
Quadrature rule for residual computation on cells.
Definition: dof_application.h:925
bool output_materials_and_levels
output_materials_and_levels if true then visualized, otherwise suppressed.
Definition: dof_application.h:493
bool print_blocks_instead_of_indices
true, if the whole blocks of FEVector solution or FEVector postprocessing to be printed to a text fil...
Definition: dof_application.h:558
Base class for applications.
Definition: application_base.h:114
unsigned int initial_refinement
Initial refinement from parameter file.
Definition: dof_application.h:890
std::vector< DataComponentInterpretation::DataComponentInterpretation > postprocessing_interpretations
postprocessing_interpretations identifies whether some postprocessing_names are scalars or parts of a...
Definition: dof_application.h:482
boost::shared_ptr< DoFHandler< dim > > dof
Pointer to the DoFHandler object.
Definition: dof_application.h:863
bool read_in_initial_solution
Bool flag used to specify if the initial solution to the problem, specially important for transient a...
Definition: dof_application.h:423
std::vector< unsigned int > solution_printing_indices
The indices of the FEVector solution to be printed to a text file.
Definition: dof_application.h:541
Point< dim > sort_direction
Direction for downstream sorting.
Definition: dof_application.h:915
IMPORTANT: Add all new solution variables and equations here !
Definition: system_management.h:271
boost::shared_ptr< FiniteElement< dim > > element
The finite element used in dof.
Definition: dof_application.h:858
Quadrature< dim-1 > quadrature_residual_bdry
Quadrature rule for residual computation on boundary faces.
Definition: dof_application.h:930
Table< 2, DoFTools::Coupling > flux_couplings
Couplings through flux bilinear forms.
Definition: dof_application.h:973
unsigned int verbosity
Controls verbosity of certain functions.
Definition: dof_application.h:991
Definition: dof_application.h:85
bool boundary_fluxes
Extend the integration loops in assemble() and residual() also to boundary faces. ...
Definition: dof_application.h:979
BlockVector< double > FEVector
The vector class used by applications.
Definition: application_data.h:39
Table< 2, DoFTools::Coupling > cell_couplings
Couplings through cell bilinear forms.
Definition: dof_application.h:967
The data type used in function calls of Application.
Definition: fe_vectors.h:59
boost::shared_ptr< Mapping< dim > > mapping
The mapping used for the transformation of mesh cells.
Definition: dof_application.h:830
double refinement_threshold
Refinement threshold for adaptive method from parameter file.
Definition: dof_application.h:895
Base class for all linear applications, i.e., all applications requiring Triangulation and DoFHandler...
Definition: dof_application.h:122
bool print_postprocessing
true, if you want to print FEVector postprocessing to a text file.
Definition: dof_application.h:532
std::vector< DataComponentInterpretation::DataComponentInterpretation > solution_interpretations
solution_interpretations identifies whether some solution_names are scalars or parts of a vector...
Definition: dof_application.h:476
std::map< unsigned int, double > boundary_constraints
List of all nodes constrained by a strong boundary condition, together with a value to be assigned...
Definition: dof_application.h:800