LCOV - code coverage report
Current view: top level - ugbase/lib_disc/spatial_disc - domain_disc.h (source / functions) Coverage Total Hit
Test: coverage.info Lines: 0.0 % 109 0
Test Date: 2025-09-21 23:31:46 Functions: 0.0 % 423 0

            Line data    Source code
       1              : /*
       2              :  * Copyright (c) 2010-2015:  G-CSC, Goethe University Frankfurt
       3              :  * Author: Andreas Vogel
       4              :  * 
       5              :  * This file is part of UG4.
       6              :  * 
       7              :  * UG4 is free software: you can redistribute it and/or modify it under the
       8              :  * terms of the GNU Lesser General Public License version 3 (as published by the
       9              :  * Free Software Foundation) with the following additional attribution
      10              :  * requirements (according to LGPL/GPL v3 §7):
      11              :  * 
      12              :  * (1) The following notice must be displayed in the Appropriate Legal Notices
      13              :  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
      14              :  * 
      15              :  * (2) The following notice must be displayed at a prominent place in the
      16              :  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
      17              :  * 
      18              :  * (3) The following bibliography is recommended for citation and must be
      19              :  * preserved in all covered files:
      20              :  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
      21              :  *   parallel geometric multigrid solver on hierarchically distributed grids.
      22              :  *   Computing and visualization in science 16, 4 (2013), 151-164"
      23              :  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
      24              :  *   flexible software system for simulating pde based models on high performance
      25              :  *   computers. Computing and visualization in science 16, 4 (2013), 165-179"
      26              :  * 
      27              :  * This program is distributed in the hope that it will be useful,
      28              :  * but WITHOUT ANY WARRANTY; without even the implied warranty of
      29              :  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
      30              :  * GNU Lesser General Public License for more details.
      31              :  */
      32              : 
      33              : #ifndef __H__UG__LIB_DISC__SPATIAL_DISC__DOMAIN_DISC__
      34              : #define __H__UG__LIB_DISC__SPATIAL_DISC__DOMAIN_DISC__
      35              : 
      36              : // other ug4 modules
      37              : #include "common/common.h"
      38              : #include "common/util/string_util.h"
      39              : 
      40              : // library intern headers
      41              : #include "subset_assemble_util.h"
      42              : #include "domain_disc_interface.h"
      43              : #include "lib_algebra/cpu_algebra_types.h"
      44              : #include "lib_disc/common/function_group.h"
      45              : #include "lib_disc/spatial_disc/elem_disc/elem_disc_assemble_util.h"
      46              : #include "lib_disc/spatial_disc/constraints/constraint_interface.h"
      47              : #include "disc_item.h"
      48              : #include "lib_disc/spatial_disc/domain_disc_interface.h"
      49              : #include "lib_disc/function_spaces/grid_function.h"
      50              : #include "lib_disc/function_spaces/error_elem_marking_strategy.h"
      51              : 
      52              : namespace ug {
      53              : 
      54              : /// \ingroup lib_disc_domain_assemble
      55              : /// @{
      56              : 
      57              : /// generic domain discretization implementing the interface
      58              : /**
      59              :  * This class template is an implementation of the IDomainDiscretization
      60              :  * interface based on the simple groupping of several local (element) 
      61              :  * discretizations and constraints.
      62              :  *
      63              :  * Functions of this class template prepare lists of the local discretizations
      64              :  * and elements (where to assemble) for every subset, whereas assembling itself
      65              :  * is performed by functions of the so-called 'global assembler class' specified
      66              :  * by the TGlobAssembler template parameter. The latter class implements
      67              :  * assembling for the generic lists of elements belonging to only one subset.
      68              :  * Cf. class StdGlobAssembler for a complete example. Note that all the functions
      69              :  * from that example should be implemented (no matter where as regular or static
      70              :  * members).
      71              :  *
      72              :  * \tparam TDomain          domain type
      73              :  * \tparam TAlgebra         algebra type
      74              :  * \tparam TGlobAssembler   global assembler type
      75              :  */
      76              : template <typename TDomain, typename TAlgebra, typename TGlobAssembler>
      77              : class DomainDiscretizationBase
      78              : :   public IDomainDiscretization<TAlgebra>,
      79              :         public IDomainMarker<TDomain>,
      80              :     protected TGlobAssembler
      81              : {
      82              :     /// Type of the global assembler
      83              :         typedef TGlobAssembler gass_type;
      84              :         
      85              :         public:
      86              :         ///     Type of Domain
      87              :                 typedef TDomain domain_type;
      88              : 
      89              :         ///     Type of algebra
      90              :                 typedef TAlgebra algebra_type;
      91              : 
      92              :         ///     Type of algebra matrix
      93              :                 typedef typename algebra_type::matrix_type matrix_type;
      94              : 
      95              :         ///     Type of algebra vector
      96              :                 typedef typename algebra_type::vector_type vector_type;
      97              : 
      98              :         /// Type of error vector
      99              :                 typedef typename CPUAlgebra::vector_type error_vector_type;
     100              : 
     101              :         ///     Type of approximation space
     102              :                 typedef ApproximationSpace<TDomain>       approx_space_type;
     103              :                 
     104              :         ///     world dimension
     105              :                 static const int dim = TDomain::dim;
     106              :                 
     107              :         public:
     108              :         ///     default Constructor
     109            0 :                 DomainDiscretizationBase(SmartPtr<approx_space_type> pApproxSpace) :
     110            0 :                         m_bErrorCalculated(false),
     111            0 :                         m_spApproxSpace(pApproxSpace), m_spAssTuner(new AssemblingTuner<TAlgebra>)
     112            0 :                 {};
     113              : 
     114              :         /// virtual destructor
     115            0 :                 virtual ~DomainDiscretizationBase() {};
     116              : 
     117              :         ///////////////////////////
     118              :         // Time independent part
     119              :         ///////////////////////////
     120              : 
     121              :         /// \copydoc IAssemble::assemble_jacobian()
     122              :                 virtual void assemble_jacobian(matrix_type& J, const vector_type& u, ConstSmartPtr<DoFDistribution> dd);
     123            0 :                 virtual void assemble_jacobian(matrix_type& J, const vector_type& u, const GridLevel& gl)
     124            0 :                 {assemble_jacobian(J, u, dd(gl));}
     125              : 
     126              :         /// \copydoc IAssemble::assemble_defect()
     127              :                 virtual void assemble_defect(vector_type& d, const vector_type& u, ConstSmartPtr<DoFDistribution> dd);
     128            0 :                 virtual void assemble_defect(vector_type& d, const vector_type& u, const GridLevel& gl)
     129            0 :                 {assemble_defect(d, u, dd(gl));}
     130              : 
     131              :         /// \copydoc IAssemble::assemble_linear()
     132              :                 virtual void assemble_linear(matrix_type& A, vector_type& b, ConstSmartPtr<DoFDistribution> dd);
     133            0 :                 virtual void assemble_linear(matrix_type& mat, vector_type& rhs, const GridLevel& gl)
     134            0 :                 {assemble_linear(mat, rhs, dd(gl));}
     135              : 
     136              :         /// \copydoc IAssemble::assemble_rhs()
     137              :                 virtual void assemble_rhs(vector_type& rhs, const vector_type& u, ConstSmartPtr<DoFDistribution> dd);
     138            0 :                 virtual void assemble_rhs(vector_type& rhs, const vector_type& u, const GridLevel& gl)
     139            0 :                 {assemble_rhs(rhs, u, dd(gl));}
     140              : 
     141              :         /// \copydoc IAssemble::assemble_rhs()
     142            0 :                 virtual void assemble_rhs(vector_type& rhs, ConstSmartPtr<DoFDistribution> dd)
     143            0 :                 {assemble_rhs(rhs, rhs, dd);}
     144            0 :                 virtual void assemble_rhs(vector_type& rhs, const GridLevel& gl)
     145            0 :                 {assemble_rhs(rhs, dd(gl));}
     146              : 
     147              :         /// \copydoc IAssemble::adjust_solution()
     148              :                 virtual void adjust_solution(vector_type& u, ConstSmartPtr<DoFDistribution> dd);
     149            0 :                 virtual void adjust_solution(vector_type& u, const GridLevel& gl)
     150            0 :                 {adjust_solution(u, dd(gl));}
     151              : 
     152              :         ///     wrapper for GridFunction
     153              :         /// \{
     154              :                 void assemble_jacobian(matrix_type& J, GridFunction<TDomain, TAlgebra>& u)
     155              :                 {assemble_jacobian(J, u, u.dof_distribution());}
     156              : 
     157              :                 void assemble_defect(vector_type& d, GridFunction<TDomain, TAlgebra>& u)
     158              :                 {assemble_defect(d, u, u.dof_distribution());}
     159              : 
     160            0 :                 void assemble_linear(matrix_type& A, GridFunction<TDomain, TAlgebra>& rhs)
     161            0 :                 {assemble_linear(A, rhs, rhs.dof_distribution());}
     162              : 
     163            0 :                 void assemble_rhs(vector_type& rhs, GridFunction<TDomain, TAlgebra>& u)
     164            0 :                 {assemble_rhs(rhs, u, u.dof_distribution());}
     165              : 
     166            0 :                 void assemble_rhs(GridFunction<TDomain, TAlgebra>& b)
     167            0 :                 {assemble_rhs(b, b.dof_distribution());}
     168              : 
     169            0 :                 void adjust_solution(GridFunction<TDomain, TAlgebra>& u)
     170            0 :                 {adjust_solution(u, u.dof_distribution());}
     171              :         /// \}
     172              : 
     173              :         ///////////////////////
     174              :         // Time dependent part
     175              :         ///////////////////////
     176              :         ///     \copydoc IDomainDiscretization::prepare_timestep()
     177              :                 virtual void prepare_timestep(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol, number future_time, ConstSmartPtr<DoFDistribution> dd);
     178            0 :                 virtual void prepare_timestep(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol, number future_time, const GridLevel& gl)
     179            0 :                 {prepare_timestep(vSol, future_time, dd(gl));}
     180              : 
     181              : 
     182              :         /// \copydoc IDomainDiscretization::prepare_timestep_elem()
     183              :                 virtual void prepare_timestep_elem(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     184              :                                               ConstSmartPtr<DoFDistribution> dd);
     185            0 :                 virtual void prepare_timestep_elem(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol, const GridLevel& gl)
     186            0 :                 {prepare_timestep_elem(vSol, dd(gl));}
     187              : 
     188              :         /// \copydoc IDomainDiscretization::assemble_jacobian()
     189              :                 virtual void assemble_jacobian(matrix_type& J,
     190              :                                                                            ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     191              :                                                                            const number s_a0,
     192              :                                                                            ConstSmartPtr<DoFDistribution> dd);
     193            0 :                 virtual void assemble_jacobian(matrix_type& J,
     194              :                                                ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     195              :                                                const number s_a, const GridLevel& gl)
     196            0 :                 {assemble_jacobian(J, vSol, s_a, dd(gl));}
     197              : 
     198              :         /// \copydoc IDomainDiscretization::assemble_defect()
     199              :                 virtual void assemble_defect(vector_type& d,
     200              :                                                                          ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     201              :                                                                          const std::vector<number>& vScaleMass,
     202              :                                                                          const std::vector<number>& vScaleStiff,
     203              :                                                                          ConstSmartPtr<DoFDistribution> dd);
     204            0 :                 virtual void assemble_defect(vector_type& d,
     205              :                                              ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     206              :                                              const std::vector<number>& vScaleMass,
     207              :                                              const std::vector<number>& vScaleStiff,
     208              :                                              const GridLevel& gl)
     209            0 :                 {assemble_defect(d, vSol, vScaleMass, vScaleStiff, dd(gl));}
     210              : 
     211              :         /// \copydoc IDomainDiscretization::assemble_linear()
     212              :                 virtual void assemble_linear(matrix_type& A, vector_type& b,
     213              :                                                                          ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     214              :                                                                          const std::vector<number>& vScaleMass,
     215              :                                                                          const std::vector<number>& vScaleStiff,
     216              :                                                                          ConstSmartPtr<DoFDistribution> dd);
     217            0 :                 virtual void assemble_linear(matrix_type& A, vector_type& b,
     218              :                                              ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     219              :                                              const std::vector<number>& vScaleMass,
     220              :                                              const std::vector<number>& vScaleStiff,
     221              :                                              const GridLevel& gl)
     222            0 :                 {assemble_linear(A, b, vSol, vScaleMass, vScaleStiff, dd(gl));}
     223              : 
     224              :         /// \copydoc IDomainDiscretization::assemble_rhs()
     225              :                 virtual void assemble_rhs(       vector_type& b,
     226              :                                                                          ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     227              :                                                                          const std::vector<number>& vScaleMass,
     228              :                                                                          const std::vector<number>& vScaleStiff,
     229              :                                                                          ConstSmartPtr<DoFDistribution> dd);
     230            0 :                 virtual void assemble_rhs(       vector_type& b,
     231              :                                              ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     232              :                                              const std::vector<number>& vScaleMass,
     233              :                                              const std::vector<number>& vScaleStiff,
     234              :                                              const GridLevel& gl)
     235            0 :                 {assemble_rhs(b, vSol, vScaleMass, vScaleStiff, dd(gl));}
     236              : 
     237              :         /// \copydoc IDomainDiscretization::adjust_solution()
     238              :                 virtual void adjust_solution(vector_type& u, number time, ConstSmartPtr<DoFDistribution> dd);
     239            0 :                 virtual void adjust_solution(vector_type& u, number time, const GridLevel& gl)
     240            0 :                 {adjust_solution(u, time, dd(gl));}
     241              : 
     242              :         ///     \copydoc IDomainDiscretization::finish_timestep()
     243              :                 virtual void finish_timestep(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol, ConstSmartPtr<DoFDistribution> dd);
     244            0 :                 virtual void finish_timestep(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol, const GridLevel& gl)
     245            0 :                 {finish_timestep(vSol, dd(gl));}
     246              : 
     247              :         /// \copydoc IDomainDiscretization::finish_timestep_elem()
     248              :                 virtual void finish_timestep_elem(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol, ConstSmartPtr<DoFDistribution> dd);
     249            0 :                 virtual void finish_timestep_elem(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol, const GridLevel& gl)
     250            0 :                 {finish_timestep_elem(vSol, dd(gl));}
     251              : 
     252              :         ///////////////////////////
     253              :         // Mass and Stiffness Matrix
     254              :         ///////////////////////////
     255              :                 virtual void init_all_exports(ConstSmartPtr<DoFDistribution> dd, bool bAsTimeDependent);
     256            0 :                 virtual void init_all_exports(const GridLevel& gl, bool bAsTimeDependent)
     257            0 :                 {init_all_exports(dd(gl), bAsTimeDependent);}
     258            0 :                 void init_all_exports(bool bAsTimeDependent)
     259            0 :                 {init_all_exports(GridLevel(GridLevel::TOP), bAsTimeDependent);}
     260              : 
     261              :         /// assembles the mass matrix
     262              :                 virtual void assemble_mass_matrix(matrix_type& M, const vector_type& u,
     263              :                                                   ConstSmartPtr<DoFDistribution> dd);
     264            0 :                 virtual void assemble_mass_matrix(matrix_type& M, const vector_type& u,
     265              :                                                   const GridLevel& gl)
     266            0 :                 {assemble_mass_matrix(M, u, dd(gl));}
     267              : 
     268              :         /// assembles the stiffness matrix
     269              :                 virtual void assemble_stiffness_matrix(matrix_type& A, const vector_type& u,
     270              :                                                        ConstSmartPtr<DoFDistribution> dd);
     271            0 :                 virtual void assemble_stiffness_matrix(matrix_type& A, const vector_type& u,
     272              :                                                        const GridLevel& gl)
     273            0 :                 {assemble_stiffness_matrix(A, u, dd(gl));}
     274              : 
     275              :         ///////////////////////////
     276              :         // Init. all exports (an optional operation, to use the exports for plotting etc.)
     277              :         ///////////////////////////
     278              : 
     279              :         ///////////////////////////////////////////////////////////
     280              :         // Error estimator                                                                              ///
     281              : public:
     282              :         /// \copydoc IDomainDiscretization::mark_error()
     283              :         // stationary
     284              :                 virtual void calc_error(const vector_type& u, ConstSmartPtr<DoFDistribution> dd, error_vector_type* u_vtk = NULL);
     285            0 :                 virtual void calc_error(const vector_type& u, const GridLevel& gl, error_vector_type* u_vtk = NULL)
     286            0 :                 {calc_error(u, dd(gl));}
     287              : 
     288            0 :                 virtual void calc_error(const GridFunction<TDomain,TAlgebra>& u)
     289            0 :                 {calc_error(u, u.dd(), NULL);}
     290            0 :                 virtual void calc_error(const GridFunction<TDomain,TAlgebra>& u, error_vector_type* u_vtk)
     291            0 :                 {calc_error(u, u.dd(), u_vtk);}
     292              : 
     293              :         // instationary
     294              :                 virtual void calc_error(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     295              :                                                                 ConstSmartPtr<DoFDistribution> dd,
     296              :                                                                 const std::vector<number>& vScaleMass,
     297              :                                                                 const std::vector<number>& vScaleStiff,
     298              :                                                                 error_vector_type* u_vtk);
     299            0 :                 virtual void calc_error(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     300              :                                                                 const std::vector<number>& vScaleMass,
     301              :                                                                 const std::vector<number>& vScaleStiff,
     302              :                                                                 const GridLevel& gl,
     303              :                                                                 error_vector_type* u_vtk)
     304              :                 {
     305            0 :                         calc_error((ConstSmartPtr<VectorTimeSeries<vector_type> >) vSol, dd(gl),
     306              :                                            vScaleMass, vScaleStiff, u_vtk);
     307            0 :                 }
     308              : 
     309              :                 virtual void mark_with_strategy(IRefiner& refiner, SmartPtr <IElementMarkingStrategy<TDomain> > strategy);
     310              : 
     311              :                 /// marks error indicators as invalid; in order to revalidate them,
     312              :                 /// they will have to be newly calculated by a call to calc_error
     313              :                 virtual void invalidate_error();
     314              : 
     315              :                 /// returns whether current error values are valid
     316              :                 virtual bool is_error_valid();
     317              : 
     318              :         protected:
     319              :                 IMultigridElementIndicators<TDomain> m_mgElemErrors;  // holding the indicators
     320              :                 bool m_bErrorCalculated;
     321              : 
     322              :         // Error estimator                                                                               //
     323              :         ///////////////////////////////////////////////////////////
     324              : 
     325              :         public:
     326              :         /// \{
     327            0 :                 virtual SmartPtr<AssemblingTuner<TAlgebra> > ass_tuner() {return m_spAssTuner;}
     328            0 :                 virtual ConstSmartPtr<AssemblingTuner<TAlgebra> > ass_tuner() const {return m_spAssTuner;}
     329              :         /// \}
     330              : 
     331              :         public:
     332              :         /// adds an element discretization to the assembling process
     333              :         /**
     334              :          * This function adds an Element Discretization to the assembling. During
     335              :          * the assembling process the elem disc will assemble a given group of
     336              :          * subsets and add the output to the passed vectors and matrices.
     337              :          * For each Elem Disc one loop over all elements of the subset will be
     338              :          * performed.
     339              :          *
     340              :          * \param[in]   elem            Element Discretization to be added
     341              :          */
     342            0 :                 void add(SmartPtr<IElemDisc<TDomain> > elem)
     343              :                 {
     344              :                 //      check that not already registered
     345            0 :                         for(size_t i = 0; i < m_vDomainElemDisc.size(); ++i)
     346            0 :                                 if(m_vDomainElemDisc[i] == elem)
     347              :                                         return;
     348              : 
     349              :                 //      set approximation space
     350            0 :                         elem->set_approximation_space(m_spApproxSpace);
     351              : 
     352              :                 //      add it
     353            0 :                         m_vDomainElemDisc.push_back(elem);
     354              :                 }
     355              : 
     356              :         /// removes a element discretization from the assembling process
     357              :         /**
     358              :          * This function removes a previously added IElemDisc from the assembling.
     359              :          * The element-wise assemblings are no longer called.
     360              :          *
     361              :          * \param[in]   elem            elem disc to be removed
     362              :          */
     363            0 :                 void remove(SmartPtr<IElemDisc<TDomain> > elem)
     364              :                 {
     365              :                         // check that already registered
     366            0 :                         for (size_t i = 0; i < m_vDomainElemDisc.size(); i++)
     367              :                         {
     368            0 :                                 if (m_vDomainElemDisc[i] == elem)
     369              :                                 {
     370              :                                         // remove constraint
     371            0 :                                         m_vDomainElemDisc.erase(m_vDomainElemDisc.begin()+i);
     372            0 :                                         return;
     373              :                                 }
     374              :                         }
     375              : 
     376              :                         UG_LOG("Tried to remove ElemDisc from DomainDisc"
     377              :                                         ", but could not find it there.");
     378              :                 }
     379              : 
     380              :                 /// adds an element error indicator to the assembling process
     381              :                 /**
     382              :                  * This function adds an Element Discretization to the assembling. During
     383              :                  * the assembling process the elem disc will assemble a given group of
     384              :                  * subsets and add the output to the passed vectors and matrices.
     385              :                  * For each Elem Disc one loop over all elements of the subset will be
     386              :                  * performed.
     387              :                  *
     388              :                  * \param[in]   elem            Element error indicator to be added
     389              :                  */
     390            0 :                 void add_elem_error_indicator(SmartPtr<IElemError<TDomain> > elem)
     391              :                 {
     392              :                         //      check that not already registered
     393            0 :                         for(size_t i = 0; i < m_vDomainElemError.size(); ++i)
     394            0 :                                 if(m_vDomainElemError[i] == elem) return;
     395              : 
     396              :                         //      set approximation space
     397            0 :                         elem->set_approximation_space(m_spApproxSpace);
     398              : 
     399              :                         //      add it
     400            0 :                         m_vDomainElemError.push_back(elem);
     401              :                 }
     402              : 
     403              :                 /// removes a element discretization from the assembling process
     404              :                 /**
     405              :                  * This function removes a previously added IElemDisc from the assembling.
     406              :                  * The element-wise assemblings are no longer called.
     407              :                  *
     408              :                  * \param[in]   elem            elem disc to be removed
     409              :                  */
     410            0 :                 void remove_elem_error_indicator(SmartPtr<IElemError<TDomain> > elem)
     411              :                 {
     412              :                         // check that 'elem' already registered
     413            0 :                         for (size_t i = 0; i < m_vDomainElemError.size(); i++)
     414              :                         {
     415            0 :                                 if (m_vDomainElemError[i] == elem)
     416              :                                 {
     417              :                                         // remove constraint
     418            0 :                                         m_vDomainElemError.erase(m_vDomainElemError.begin()+i);
     419            0 :                                         return;
     420              :                                 }
     421              :                         }
     422              : 
     423              :                         UG_LOG("Tried to remove ElemError from DomainDisc"
     424              :                                         ", but could not find it there.");
     425              :                 }
     426              : 
     427              : 
     428              :         /// adds a constraint to the assembling process
     429              :         /**
     430              :          * This function adds a IConstraint to the assembling. The constraint is
     431              :          * called when all element-wise assembling have been performed.
     432              :          *
     433              :          * \param[in]   pp              Constraint to be added
     434              :          */
     435            0 :                 void add(SmartPtr<IDomainConstraint<TDomain, TAlgebra> > pp)
     436              :                 {
     437              :                 //      check that not already registered
     438            0 :                         for(size_t i = 0; i < m_vConstraint.size(); ++i)
     439            0 :                                 if(m_vConstraint[i] == pp)
     440              :                                         return;
     441              : 
     442              :                 //      set approximation space
     443            0 :                         pp->set_approximation_space(m_spApproxSpace);
     444              : 
     445              :                 //      add constraint
     446            0 :                         m_vConstraint.push_back(pp);
     447              :                 }
     448              : 
     449              :         /// removes a constraint from the assembling process
     450              :         /**
     451              :          * This function removes a previously added IConstraint from the assembling.
     452              :          * The constraint's assemblings are no longer called.
     453              :          *
     454              :          * \param[in]   pp              constraint to be removed
     455              :          */
     456            0 :                 void remove(SmartPtr<IDomainConstraint<TDomain, TAlgebra> > pp)
     457              :                 {
     458              :                         // check that already registered
     459            0 :                         for (size_t i = 0; i < m_vConstraint.size(); i++)
     460              :                         {
     461            0 :                                 if (m_vConstraint[i] == pp)
     462              :                                 {
     463              :                                         // remove constraint
     464            0 :                                         m_vConstraint.erase(m_vConstraint.begin()+i);
     465            0 :                                         return;
     466              :                                 }
     467              :                         }
     468              : 
     469              :                         UG_LOG("Tried to remove DomainConstraint from DomainDisc"
     470              :                                         ", but could not find it there.");
     471              :                 }
     472              : 
     473              :         /// adds a disc item to the assembling process
     474              :         /**
     475              :          * This function adds a IDiscretizationItem to the assembling. The contained
     476              :          * elem discs and constraints are sorted into the lists
     477              :          *
     478              :          * \param[in]   di              Disc Item
     479              :          */
     480            0 :                 void add(SmartPtr<IDiscretizationItem<TDomain, TAlgebra> > di)
     481              :                 {
     482              :                 //      add elem discs
     483            0 :                         for(size_t i = 0; i < di->num_elem_disc(); ++i)
     484            0 :                                 add(di->elem_disc(i));
     485              : 
     486              :                 //      add constraints
     487            0 :                         for(size_t i = 0; i < di->num_constraint(); ++i)
     488            0 :                                 add(di->constraint(i));
     489            0 :                 }
     490              : 
     491              :         ///     returns number of registered constraints
     492            0 :                 virtual size_t num_constraints() const {return m_vConstraint.size();}
     493              : 
     494              :         ///     returns the i'th constraint
     495            0 :                 virtual SmartPtr<IConstraint<TAlgebra> > constraint(size_t i) {return m_vConstraint[i];}
     496              : 
     497            0 :                 SmartPtr<approx_space_type> approximation_space ()                                {return m_spApproxSpace;}
     498            0 :                 ConstSmartPtr<approx_space_type> approximation_space () const     {return m_spApproxSpace;}
     499              : 
     500              :         protected:
     501              :         ///     set the approximation space in the elem discs and extract IElemDiscs
     502              :                 void update_elem_discs();
     503              :                 void update_elem_errors();
     504              :                 void update_constraints();
     505              :                 void update_disc_items();
     506              :                 void update_error_items();
     507              : 
     508              :         protected:
     509              :         ///     returns the level dof distribution
     510            0 :                 ConstSmartPtr<DoFDistribution> dd(const GridLevel& gl) const{return m_spApproxSpace->dof_distribution(gl);}
     511              : 
     512              :         protected:
     513              :         ///     vector holding all registered elem discs
     514              :                 std::vector<SmartPtr<IElemDisc<TDomain> > > m_vDomainElemDisc;
     515              : 
     516              :         ///     vector holding all registered elem discs
     517              :                 std::vector<IElemDisc<TDomain>*> m_vElemDisc;
     518              : 
     519              :         ///     vector holding all registered elem discs
     520              :                 std::vector<SmartPtr<IElemError<TDomain> > > m_vDomainElemError;
     521              : 
     522              :         ///     vector holding all registered elem discs
     523              :                 std::vector<IElemError<TDomain>*> m_vElemError;
     524              : 
     525              :         ///     vector holding all registered constraints
     526              :                 std::vector<SmartPtr<IDomainConstraint<TDomain, TAlgebra> > > m_vConstraint;
     527              : 
     528              :         ///     current approximation space
     529              :                 SmartPtr<approx_space_type> m_spApproxSpace;
     530              :                 
     531              :         ///     this object provides tools to adapt the assemble routine
     532              :                 SmartPtr<AssemblingTuner<TAlgebra> > m_spAssTuner;
     533              :         
     534              :         private:
     535              :         //---- Auxiliary function templates for the assembling ----//
     536              :         //      These functions call the corresponding functions from the global assembler for a composed list of elements:
     537              :         //-- for stationary problems --//
     538              :         template <typename TElem>
     539              :         void AssembleMassMatrix(                const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     540              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     541              :                                                                         int si, bool bNonRegularGrid,
     542              :                                                                         matrix_type& M,
     543              :                                                                         const vector_type& u);
     544              :         template <typename TElem>
     545              :         void AssembleStiffnessMatrix(   const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     546              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     547              :                                                                         int si, bool bNonRegularGrid,
     548              :                                                                         matrix_type& A,
     549              :                                                                         const vector_type& u);
     550              :         template <typename TElem>
     551              :         void AssembleJacobian(                  const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     552              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     553              :                                                                         int si, bool bNonRegularGrid,
     554              :                                                                         matrix_type& J,
     555              :                                                                         const vector_type& u);
     556              :         template <typename TElem>
     557              :         void AssembleDefect(                    const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     558              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     559              :                                                                         int si, bool bNonRegularGrid,
     560              :                                                                         vector_type& d,
     561              :                                                                         const vector_type& u);
     562              :         template <typename TElem>
     563              :         void AssembleLinear(                    const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     564              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     565              :                                                                         int si, bool bNonRegularGrid,
     566              :                                                                         matrix_type& A,
     567              :                                                                         vector_type& rhs);
     568              :         template <typename TElem>
     569              :         void AssembleRhs(                               const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     570              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     571              :                                                                         int si, bool bNonRegularGrid,
     572              :                                                                         vector_type& rhs,
     573              :                                                                         const vector_type& u);
     574              :         template <typename TElem>
     575              :         void AssembleErrorEstimator(    const std::vector<IElemError<domain_type>*>& vElemDisc,
     576              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     577              :                                                                         int si, bool bNonRegularGrid,
     578              :                                                                         const vector_type& u);
     579              :         //-- for time-dependent problems --//
     580              :         template <typename TElem>
     581              :         void PrepareTimestepElem(               const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     582              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     583              :                                                                         int si, bool bNonRegularGrid,
     584              :                                                                         ConstSmartPtr<VectorTimeSeries<vector_type> > vSol);
     585              :         template <typename TElem>
     586              :         void AssembleJacobian(                  const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     587              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     588              :                                                                         int si, bool bNonRegularGrid,
     589              :                                                                         matrix_type& J,
     590              :                                                                         ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     591              :                                                                         number s_a0);
     592              :         template <typename TElem>
     593              :         void AssembleDefect(                    const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     594              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     595              :                                                                         int si, bool bNonRegularGrid,
     596              :                                                                         vector_type& d,
     597              :                                                                         ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     598              :                                                                         const std::vector<number>& vScaleMass,
     599              :                                                                         const std::vector<number>& vScaleStiff);
     600              :         template <typename TElem>
     601              :         void AssembleLinear(                    const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     602              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     603              :                                                                         int si, bool bNonRegularGrid,
     604              :                                                                         matrix_type& A,
     605              :                                                                         vector_type& rhs,
     606              :                                                                         ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     607              :                                                                         const std::vector<number>& vScaleMass,
     608              :                                                                         const std::vector<number>& vScaleStiff);
     609              :         template <typename TElem>
     610              :         void AssembleRhs(                               const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     611              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     612              :                                                                         int si, bool bNonRegularGrid,
     613              :                                                                         vector_type& rhs,
     614              :                                                                         ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
     615              :                                                                         const std::vector<number>& vScaleMass,
     616              :                                                                         const std::vector<number>& vScaleStiff);
     617              :         template <typename TElem>
     618              :         void FinishTimestepElem(                        const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     619              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     620              :                                                                         int si, bool bNonRegularGrid,
     621              :                                                                         ConstSmartPtr<VectorTimeSeries<vector_type> > vSol);
     622              :         template <typename TElem>
     623              :         void AssembleErrorEstimator(    const std::vector<IElemError<domain_type>*>& vElemDisc,
     624              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     625              :                                                                         int si, bool bNonRegularGrid,
     626              :                                                                         const std::vector<number>& vScaleMass,
     627              :                                                                         const std::vector<number>& vScaleStiff,
     628              :                                                                         ConstSmartPtr<VectorTimeSeries<vector_type> > vSol);
     629              :         template <typename TElem>
     630              :         void InitAllExports(                    const std::vector<IElemDisc<domain_type>*>& vElemDisc,
     631              :                                                                         ConstSmartPtr<DoFDistribution> dd,
     632              :                                                                         int si, bool bNonRegularGrid, bool bAsTimeDependent);
     633              : };
     634              : 
     635              : /// domain discretization implementing the interface
     636              : /**
     637              :  * This class template is an implementation of the IDomainDiscretization
     638              :  * interface based on the simple groupping of several local (element) 
     639              :  * discretizations and constraints.
     640              :  *
     641              :  * This is the "usual" and "trivial" global discretizations: all the
     642              :  * local (element) discretizations are perfomed once for every element, and
     643              :  * their contributions are algebraically added to the global data.
     644              :  *
     645              :  * \tparam TDomain          domain type
     646              :  * \tparam TAlgebra         algebra type
     647              :  */
     648              : template <typename TDomain, typename TAlgebra>
     649              : class DomainDiscretization
     650              : :       public DomainDiscretizationBase<TDomain, TAlgebra, StdGlobAssembler<TDomain, TAlgebra> >
     651              : {
     652              :     /// Type of the global assembler
     653              :         typedef StdGlobAssembler<TDomain, TAlgebra> gass_type;
     654              :         
     655              :         public:
     656              :         ///     Type of Domain
     657              :                 typedef TDomain domain_type;
     658              : 
     659              :         ///     Type of algebra
     660              :                 typedef TAlgebra algebra_type;
     661              : 
     662              :         ///     Type of algebra matrix
     663              :                 typedef typename algebra_type::matrix_type matrix_type;
     664              : 
     665              :         ///     Type of algebra vector
     666              :                 typedef typename algebra_type::vector_type vector_type;
     667              : 
     668              :         ///     Type of approximation space
     669              :                 typedef ApproximationSpace<TDomain>       approx_space_type;
     670              :                 
     671              :         ///     world dimension
     672              :                 static const int dim = TDomain::dim;
     673              :                 
     674              :         public:
     675              :         ///     default Constructor
     676            0 :                 DomainDiscretization(SmartPtr<approx_space_type> pApproxSpace)
     677            0 :                 : DomainDiscretizationBase<domain_type, algebra_type, gass_type> (pApproxSpace)
     678            0 :                 {};
     679              : 
     680              :         /// virtual destructor
     681            0 :                 virtual ~DomainDiscretization() {};
     682              : };
     683              : 
     684              : /// @}
     685              : 
     686              : } // end namespace ug
     687              : 
     688              : // include documentation
     689              : #include "domain_disc_impl.h"
     690              : 
     691              : #endif /* __H__UG__LIB_DISC__SPATIAL_DISC__DOMAIN_DISC__ */
        

Generated by: LCOV version 2.0-1