OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
matrix_application.h
Go to the documentation of this file.
1 //---------------------------------------------------------------------------
2 // $Id: matrix_application.h 1348 2013-08-16 00:45:07Z secanell $
3 //
4 // Copyright (C) 2006, 2007, 2008, 2009 by Guido Kanschat
5 //
6 // This file is subject to QPL and may not be distributed
7 // without copyright and license information. Please refer
8 // to the file deal.II/doc/license.html for the text and
9 // further information on this license.
10 //
11 //---------------------------------------------------------------------------
12 
13 #ifndef __deal2__appframe__matrix_application_h
14 #define __deal2__appframe__matrix_application_h
15 
16 #include <base/quadrature.h>
17 #include <base/thread_management.h>
18 #include <base/multithread_info.h>
19 #include <lac/full_matrix.h>
20 #include <lac/sparse_matrix.h>
21 #include <lac/block_sparsity_pattern.h>
22 #include <grid/tria.h>
23 #include <grid/tria_iterator.h>
24 #include <grid/tria_accessor.h>
25 #include <dofs/dof_tools.h>
26 #include <base/mg_level_object.h>
28 #include <appframe/matrix_base.h>
29 #include <boost/shared_ptr.hpp>
30 
31 namespace AppFrame
32 {
103  template <int dim>
105  public DoFApplication<dim>
106  {
107  public:
117 
123 
129 
131 
132 
143  MatrixApplication(boost::shared_ptr<ApplicationData> data =
144  boost::shared_ptr<ApplicationData>(),
145  bool multigrid = false);
146 
152  MatrixApplication(bool multigrid);
153 
166  bool triangulation_only,
167  bool multigrid = false);
168 
181  void _initialize(ParameterHandler& param);
182 
193  virtual void declare_parameters(ParameterHandler& param);
194 
200  virtual void initialize (ParameterHandler& param);
202 
203 
204 
218  void add_matrix(unsigned int row, unsigned int col);
219 
225  void add_mg_matrix(unsigned int row, unsigned int col);
226 
234  void add_vector_for_cell_matrix(std::string name,
235  unsigned int block,
236  unsigned int components,
237  bool values,
238  bool derivatives);
247  void add_vector_for_flux_matrix(std::string name,
248  unsigned int block,
249  unsigned int components,
250  bool values,
251  bool derivatives);
253 
261  virtual void post_cell_assemble();
262 
264 
265 
271  void remesh_matrices();
272 
279  virtual void remesh();
280 
290  void assemble(const FEVectors&);
292 
293 
294 
301  virtual void cell_matrix(MatrixVector& cell_matrices,
302  const typename DoFApplication<dim>::CellInfo& cell);
303 
308  virtual void bdry_matrix(MatrixVector& face_matrices,
309  const typename DoFApplication<dim>::FaceInfo& face);
310 
315  virtual void face_matrix(MatrixVector& matrices11,
316  MatrixVector& matrices12,
317  MatrixVector& matrices21,
318  MatrixVector& matrices22,
319  const typename DoFApplication<dim>::FaceInfo& face1,
320  const typename DoFApplication<dim>::FaceInfo& face2);
321 
330  virtual void mg_cell_matrix(MatrixVector& cell_matrices,
331  const typename DoFApplication<dim>::CellInfo& cell);
332 
339  virtual void mg_bdry_matrix(MatrixVector& face_matrices,
340  const typename DoFApplication<dim>::FaceInfo& face);
341 
348  virtual void mg_face_matrix(MatrixVector& matrices11,
349  MatrixVector& matrices12,
350  MatrixVector& matrices21,
351  MatrixVector& matrices22,
352  const typename DoFApplication<dim>::FaceInfo& face1,
353  const typename DoFApplication<dim>::FaceInfo& face2);
355  template <class ASSEMBLER>
356  class LocalMatrixIntegrator : public ASSEMBLER,
358  {
359  SmartPointer<MatrixApplication<dim> > app;
360  bool multigrid;
361  public:
362  LocalMatrixIntegrator(bool mg = false)
363  : multigrid(mg)
364  {}
365 
367  const MeshWorker::BlockInfo& bi, std::vector<typename ASSEMBLER::MatrixPtr>& m)
368  {
369  app = a;
370  ASSEMBLER::initialize(bi, m);
371  }
372 
377  {
378  if (multigrid)
379  app->mg_cell_matrix(this->M11, cell);
380  else
381  app->cell_matrix(this->M11, cell);
382  assemble(cell);
383  }
384 
385 
390  {
391  if (multigrid)
392  app->mg_bdry_matrix(this->M11, face);
393  else
394  app->bdry_matrix(this->M11, face);
395  assemble(face);
396  }
397 
401  void face(const typename DoFApplication<dim>::FaceInfo& face1,
402  const typename DoFApplication<dim>::FaceInfo& face2)
403  {
404  if (multigrid)
405  app->mg_face_matrix(this->M11, this->M12, this->M21, this->M22, face1, face2);
406  else
407  app->face_matrix(this->M11, this->M12, this->M21, this->M22, face1, face2);
408  assemble(face1, face2);
409  }
410  };
411 
412  protected:
441  Table<2, DoFTools::Coupling> cell_couplings;
442 
447  Table<2, DoFTools::Coupling> flux_couplings;
448 
453  Table<2, DoFTools::Coupling> mg_cell_couplings;
454 
459  Table<2, DoFTools::Coupling> mg_flux_couplings;
460 
465  boost::shared_ptr<Quadrature<dim> > quadrature_assemble_cell;
466 
471  boost::shared_ptr<Quadrature<dim-1> > quadrature_assemble_face;
472 
473  private:
477  BlockSparsityPattern sparsity;
481  MGLevelObject<BlockSparsityPattern> mg_sparsity;
482 
487  MGLevelObject<BlockSparsityPattern> mg_sparsity_interface;
488 
489  protected:
503  std::vector<boost::shared_ptr<MatrixType> > matrices;
504 
508  std::vector<boost::shared_ptr<MGMatrixType> > mg_matrices;
509 
516  std::vector<boost::shared_ptr<MGMatrixType> > mg_interface_matrices_up;
517 
524  std::vector<boost::shared_ptr<MGMatrixType> > mg_interface_matrices_down;
525 
526  private:
531 
539  };
540 
541 }
542 
543 #endif